Full API

Autogenerated API list

WaveguideQED.LazySumKetType
LazySumKet(kets...)

Lazy sum of LazyTensorKets that is used to express entanglement between subsystems in LazyTensorKets.

source
WaveguideQED.OnePhotonViewType
OnePhotonView{T} <: AbstractVector{T}

Structure for viewing onephoton excitations in waveguides of the form $W^\dagger(\xi) |0 \rangle = \int_{t_0}^{t_{end}} dt \xi(t) w_{\mathrm{i}}^\dagger(t) |\emptyset \rangle$ where $i$ is the index of the waveguide. See onephoton on how to create onephoton wavepackets and view_waveguide on how to index when there are multiple systems.

Examples

using QuantumOptics;
times = 0:1:10;
bw = WaveguideBasis(1,times);
ψ1 = onephoton(bw,x->1,norm=false);
OnePhotonView(ψ1) == ones(length(times))
bc = FockBasis(2);
ψ1Cavity = fockstate(bc,2,norm=false) ⊗ ψ1;
OnePhotonView(ψCavity,[3,:]) == ones(length(times))
bw =  WaveguideBasis(1,3,times);
ψ2 = onephoton(bw,2,x->1,norm=false)
OnePhotonView(ψ,2) == ones(length(times))
ψ2Cavity = fockstate(bc,2) ⊗ ψ2;
OnePhotonView(ψ2Cavity,2,[3,:]) == ones(length(times))
source
WaveguideQED.OnePhotonViewMethod
OnePhotonView(ψ::T) where {T<:SingleWaveguideKet}

ψ contains only a single waveguide and no waveguide index is needed. No index of the system is provided and groundstate is assumed. Thus returns OnePhotonView(ψ,index) with index = [1,:,1,...] with : at the location of the waveguide and 1 in every other position.

source
WaveguideQED.OnePhotonViewMethod
OnePhotonView(ψ::T,WI::Int,index::I) where {T<:MultipleWaveguideKet,I<:Union{Vector{Any},Vector{Int64},Tuple{Vararg{Union{Int64,Colon}}}}}

ψ contains only a multiple waveguides and waveguide index WI is needed. index should follow syntax outlined in view_waveguide.

source
WaveguideQED.OnePhotonViewMethod
OnePhotonView(ψ::T,index::I) where {T<:SingleWaveguideKet,I<:Union{Vector{Any},Vector{Int64},Tuple{Vararg{Union{Int64,Colon}}}}}

ψ contains only a single waveguide and no waveguide index is needed. index should follow syntax outlined in view_waveguide.

source
WaveguideQED.OnePhotonViewMethod
OnePhotonView(ψ::T,WI::Int)  where {T<:MultipleWaveguideKet}

ψ contains only a multiple waveguides and waveguide index WI is needed. No index of the system is provided and groundstate is assumed (see view_waveguide).

source
WaveguideQED.TwoPhotonViewType
TwophotonView{T} <: AbstractMatrix{T}

Structure for viewing twophoton state of the form $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{j}}^\dagger(t') |\emptyset \rangle$ where $i$ and $j$ are the indeces of the waveguide as matrix. See also twophoton and view_waveguide.

Examples

Basic viewing:

using LinearAlgebra; #hide
times = 1:1:10;
bw = WaveguideBasis(2,times);
ψ = twophoton(bw,(t1,t2)->1,norm=false);
ψview = TwoPhotonView(ψ);
ψview == ones((length(times),length(times)))

Viewing state combined with other basis:

using QuantumOptics;
bc = FockBasis(2);
ψcombined = fockstate(bc,2) ⊗ ψ;
ψview = TwoPhotonView(ψcombined,[3,:]);
ψview == ones((length(times),length(times)))

Viewing twophoton state in waveguide 2 with multiple waveguides

bw = WaveguideBasis(2,3,times)
ψ = twophoton(bw,2,(t1,t2)->1,norm=false);
ψview = TwoPhotonView(ψ,2);
ψview == ones((length(times),length(times)))

Viewing twophoton state in waveguide 2 with multiple waveguides combined with other basis:

ψcombined = fockstate(bc,2) ⊗ ψ;
ψview = TwoPhotonView(ψcombined,2,[3,:]);
ψview == ones((length(times),length(times)))

Viewing twophotons across waveguide 1 and 2

bw = WaveguideBasis(2,3,times)
ψ = twophoton(bw,1,2,(t1,t2)->1,norm=false);
ψview = TwoPhotonView(ψ,1,2);
ψview == ones((length(times),length(times)))

Viewing twophotons across waveguide 1 and 2 combined with other basis:

ψcombined = fockstate(bc,2) ⊗ ψ;
ψview = TwoPhotonView(ψcombined,1,2,[3,:]);
ψview == ones((length(times),length(times)))
source
WaveguideQED.TwoPhotonViewMethod
TwoPhotonView(ψ::T,WI::F,index::I)  where {T<:MultipleWaveguideKet,I<:Union{Vector{Any},Vector{Int64},Tuple{Vararg{Union{Int64,Colon}}}},F<:Union{Vector{Int64},Tuple{Vararg{Int64}}}}

State ψ contains multiple waveguides.

Waveguide indeces provided as tuple or vector of length 2, means viewing viewing the state $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{i}}^\dagger(t') |\emptyset \rangle$ with i = WI[1] and j = WI[2].

Index should follow syntax highlighted in view_waveguide.

source
WaveguideQED.TwoPhotonViewMethod
TwoPhotonView(ψ::T,WI1::Int,WI2::Int,index::I)  where {T<:MultipleWaveguideKet,I<:Union{Vector{Any},Vector{Int64},Tuple{Vararg{Union{Int64,Colon}}}}}

State ψ contains multiple waveguides.

Two waveguide indeces means viewing viewing the state $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{i}}^\dagger(t') |\emptyset \rangle$ with i = WI1 and j = WI2.

Index should follow syntax highlighted in view_waveguide.

source
WaveguideQED.TwoPhotonViewMethod
TwoPhotonView(ψ::T,WI::Int,index::I)  where {T<:MultipleWaveguideKet,I<:Union{Vector{Any},Vector{Int64},Tuple{Vararg{Union{Int64,Colon}}}}}

State ψ contains multiple waveguides.

Only one waveguide index i=WI means viewing the state $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{i}}^\dagger(t') |\emptyset \rangle$.

Index should follow syntax highlighted in view_waveguide.

source
WaveguideQED.TwoPhotonViewMethod
TwoPhotonView(ψ::T,index::I) where {T <: SingleWaveguideKet,I<:Union{Vector{Any},Vector{Int64},Tuple{Vararg{Union{Int64,Colon}}}}}

State ψ only contains one waveguide. Index should follow syntax highlighted in view_waveguide.

source
WaveguideQED.TwoPhotonViewMethod
TwoPhotonView(ψ::T,WI::I) where {T <: MultipleWaveguideKet,I<:Union{Vector{Int64},Tuple{Vararg{Int64}}}}

State ψ contains multiple waveguides.

Waveguide indeces provided as tuple or vector of length 2, means viewing viewing the state $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{i}}^\dagger(t') |\emptyset \rangle$ with i = WI[1] and j = WI[2].

No index provided so groundstate is assumed. See view_waveguide.

source
WaveguideQED.TwoPhotonViewMethod
TwoPhotonView(ψ::T,WI1::Int,WI2::Int) where {T <: MultipleWaveguideKet}

State ψ contains multiple waveguides.

Two waveguide indeces means viewing viewing the state $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{i}}^\dagger(t') |\emptyset \rangle$ with i = WI1 and j = WI2.

No index provided so groundstate is assumed. See view_waveguide.

source
WaveguideQED.TwoPhotonViewMethod
TwoPhotonView(ψ::T,WI::Int) where {T <: MultipleWaveguideKet}

State ψ contains multiple waveguides and waveguide index WI required.

Only one waveguide index i=WI means viewing the state $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{i}}^\dagger(t') |\emptyset \rangle$.

No index provided so groundstate is assumed. See view_waveguide.

source
WaveguideQED.WaveguideBasisType
WaveguideBasis(Np,Nw, times)
WaveguideBasis(Np, times)

Basis for time binned Waveguide where Np is the number of photons in the waveguide and Nw the number of waveguides (default is 1). Currently number of photons is restricted to either 1 or 2. Times is timeinterval over which the photon state should be binned.

source
WaveguideQED.WaveguideCreateType
WaveguideCreate{B1,B2,N,idx} <: WaveguideOperator{B1,B2}

Operator structure for dispatching creation operation on Waveguide state. Np is used to dispatch one or two photon routine and idx denotes the index of the waveguide the operator is acting on.

source
WaveguideQED.WaveguideDestroyType
WaveguideDestroy{B1,B2,Np,idx} <: WaveguideOperator{B1,B2}

Operator structure for dispatching annihilation operation on Waveguide state. Np is used to dispatch one or two photon routine and idx denotes the index of the waveguide the operator is acting on.

source
WaveguideQED.WaveguideTransformType
WaveguideTransform{B1,B2,Np,idx} <: WaveguideOperator{B1,B2}

Operator structure for transforming an output state $\ket{\psi}_{\mathrm{out}} = \mathbf{C} \ket{\tilde{\psi}}_{\mathrm{out}}$ Np is used to dispatch one or two photon routine and idx denotes the index of the waveguide the operator is acting on. See also effective_hamiltonian.

Examples:

times = 0:0.1:10
bw = WaveguideBasis(1,2,times)
C = 1/sqrt(2)*[1.0 -im;-im 1]
C_transform = WaveguideTransform(bw,C)
psi_tilde = ket(bw)
psi = C_transform*psi_tilde
source
QuantumInterface.tensorMethod
tensor(a::AbstractOperator,b::CavityWaveguideAbsorption)
tensor(a::CavityWaveguideAbsorption,b::AbstractOperator)
tensor(a::AbstractOperator,b::CavityWaveguideEmission)
tensor(a::CavityWaveguideEmission,b::AbstractOperator)

Methods for tensorproducts between QuantumOptics.jl operators and CavityWaveguideOperator.

source
QuantumOpticsBase.createMethod
create(basis::WaveguideBasis{Np,1}) where {Np}
create(basis::WaveguideBasis{Np,Nw},i::Int) where {Np,Nw}

Creation operator $w^\dagger$ for WaveguideBasis $w_i(t_k)^\dagger | \emptyset \rangle = | 1_k \emptyset \rangle_i$ with cutoff (maximum number of photons Np) where i is the index of the waveguide. $t_k$ is determined by the timeindex property of the operator which can be changed by set_waveguidetimeindex!(a::WaveguideOperator,k::Int)

Arguments

  • basis of type WaveguideBasis, defines cutoff photon number Np and number of waveguides Nw
  • i determines which waveguide the operator acts on and should be i ≤ Nw. If Nw=1 then i=1 is assumed (there is only one waveguide).

Returns

WaveguideCreate

source
QuantumOpticsBase.destroyMethod
destroy(basis::WaveguideBasis{Np,1}) where {Np}
destroy(basis::WaveguideBasis{Np,Nw},i::Int) where {Np,Nw}

Annihilation operator $w$ for WaveguideBasis $w_i(t_k) | 1_j \emptyset \rangle_i = \delta_{k,j} | \emptyset \rangle$ with cutoff (maximum number of photons Np) where i is the index of the waveguide. $t_k$ is determined by the timeindex property of the operator which can be changed by set_waveguidetimeindex!(a::WaveguideOperator,k::Int)

Arguments

  • basis of type WaveguideBasis, defines cutoff photon number Np and number of waveguides Nw
  • i determines which waveguide the operator acts on and should be i ≤ Nw. If Nw=1 then i=1 is assumed (there is only one waveguide).

Returns

WaveguideDestroy

source
WaveguideQED.detect_double_clickMethod
detect_double_click(ψ,detector1,detector2,projection)
detect_double_click(ψ,detector1,detector2)

Calculate probability of observing projection after beamsplitter operation and two subsequent detection events defined in detector1 and detector2 on the state ψ.

Arguments

  • ψ can be either LazyTensorKet or LazySumKet and is the state on which the beamsplitter and detection is applied
  • detector1 defines the first beamsplitter and subsequent detection operation. See Detector for more details on how to define.
  • detector2 defines the second beamsplitter and subsequent detection operation. See Detector for more details on how to define.
  • projection if given is a LazyTensorKet or LazySumKet which projects onto the state after the measurement. If no projection is given, instead the total probability of having the detector click is given by applying all possible combinations of projections using get_all_projectors.

Returns

  • If projection is given: Returns probability of having detector1 and detector2 click and being in state defined by projection
  • If projection is not given: Returns the total probability of having detector1 and detector2 click by applying all possibile projections with zerophotons in the waveguide using get_all_projectors.

See Beamsplitter for an example on how to use.

source
WaveguideQED.detect_single_clickMethod
detect_single_click(ψ,detector::Detector,projection)
detect_single_click(ψ,detector::Detector)

Calculate probability of observing projection after beamsplitter operation and subsequent detection event defined in detector on the state ψ.

Arguments

  • ψ can be either LazyTensorKet or LazySumKet and is the state on which the beamsplitter and detection is applied
  • detector defines the beamsplitter and subsequent detection operation. See Detector for more details on how to define.
  • projection if given is a LazyTensorKet or LazySumKet which projects onto the state after the measurement. If no projection is given, instead the total probability of having the detector click is given by applying all possible combinations of projections using get_all_projectors.

Returns

  • If projection is given returns probability of having detector click and being in state defined by projection
  • If projection is not given returns the total probability of having a the detector click (only a single click, for double clicks use detect_double_click) by applying all possibile projections with zerophotons in the waveguide using get_all_projectors.

See Beamsplitter for an example on how to use.

source
WaveguideQED.effective_hamiltonianMethod
effective_hamiltonian(bw::WaveguideBasis{Np,Nw},bs,C,G) where {Np,Nw}

Return the effective Hamiltonian for a waveguide system with basis bw coupled to a local emitter/cavity with the basis bc. C and G here determine the input output relations of the coupling such that: $\frac{d}{d t} a=-i\left[a, H_{\mathrm{c}}\right]-\Sigma a+\mathbf{k}^T \mathbf{W}_{\mathrm{in}}$, and $\mathbf{W}_{\text {out }}(t)=\mathbf{C} \mathbf{W}_{\mathrm{in}}(t)+a(t)\mathbf{d}$. Here $\tilde{\mathbf{d}} = \tilde{\mathbf{k}} = -i \left(\left(\mathbf{I}+\frac{i}{2} \mathbf{V}\right)^{-1} \right) \mathbf{\Gamma}$, where $\mathbf{\Gamma}$ is determined by the vector G with a length equal to the number of waveguide in bw. $\mathbf{C}$ is here given by C which is of dimensions (nw,nw) with nw being the number of waveguides in bw.

The resulting state $\ket{\tilde{\psi}}_{\mathrm{out}}$ from simulating with this Hamiltonian needs to be transformed using WaveguideTransform such that $\ket{\psi}_{\mathrm{out}} = \mathbf{C} \ket{\tilde{\psi}}_{\mathrm{out}}$ to get the correct input-output relations.

Examples:

times = 0:0.1:10
dt = times[2] - times[1]

bw = WaveguideBasis(1,2,times)
be = FockBasis(1)


C = [0 -1.0im; -1.0im 0]
γ=1
G = [sqrt(γ/dt),sqrt(γ/dt)]

H = effective_hamiltonian(bw,be,C,G)

ξfun(t1,σ,t0) = sqrt(2/σ)* (log(2)/pi)^(1/4)*exp(-2*log(2)*(t1-t0)^2/σ^2)

ψ_in = onephoton(bw,1,ξfun,1,5) ⊗ fockstate(be,0)
ψ_out_tilde = waveguide_evolution(times,ψ_in,H)


C_transform = WaveguideTransform(bw,C) ⊗ identityoperator(be)
ψ_out = C_transform*ψ_out_tilde
source
WaveguideQED.expect_waveguideMethod
expect_waveguide(O, psi)

Returns the sum of expectation values of O over waveguide time indeces. This is usefull to extract information from waveguide states. For example the number of photons can be computed as: expect_waveguide(create(bw)*destroy(bw), psi, times), where BW is the waveguidebasis. Mathematically this is equivalent to $sum_n \langle \psi | O(t_n) | \psi \rangle$. For $O(t_k) = w_k^\dagger w_k$ acting on a single photon state $|\psi \rangle = \sum_k \xi(t_k) w_k^\dagger |0\rangle$, one would therefore have: $\langle O \rangle = \sum_k |\xi(t_k)|^2 = 1$, where the last follows from $\xi(t)$ being normalized.

source
WaveguideQED.fftketMethod
fftket(psi::Ket,idx)

Return the FFT $\xi(\omega)$ of a onephoton wavefunction $\xi(t)$. Idx determines which waveguide.

source
WaveguideQED.get_waveguide_locationMethod
get_waveguide_location(basis::WaveguideBasis)
get_waveguide_location(basis::CompositeBasis)

Return index of WaveguideBasis location in Hilbert space of basis b. Btotal = waveguidebasis ⊗ otherbasis where waveguidebasis is a WaveguideBasis and otherbasis some other basis then get_waveguide_location(Btotal) returns 1. While Btotal = otherbasis ⊗ waveguidebasis with get_waveguide_location(Btotal) returns 2.

source
WaveguideQED.get_waveguide_operatorsMethod
get_waveguide_operators(basis::LazySum)
get_waveguide_operators(basis::LazyProduct)
get_waveguide_operators(basis::LazyTensor)
get_waveguide_operators(basis::Tuple)
get_waveguide_operators(basis::Array)
get_waveguide_operators(basis::WaveguideOperator)

Returns all WaveguideOperator in LazyOperator or from a list of operators. If no WaveguideOperator is found, and empty array is returned.

source
WaveguideQED.onephotonFunction

Create a onephoton wavepacket in waveguide of the form $W^\dagger(\xi) |0 \rangle = \int_{t_0}^{t_{end}} dt \xi(t) w_{\mathrm{i}}^\dagger(t) |\emptyset \rangle = \sum_{k}$ where $i$ is the index of the waveguide and return it as a Ket. See also WaveguideBasis and OnePhotonView for how to view the state.

Examples

dt = 1
times = 1:dt:10;
bw = WaveguideBasis(1,times);
ψ = onephoton(bw,x->dt,norm=false);
OnePhotonView(ψ) == ones(length(times))
vec = collect(1:1:10);
ψ = onephoton(bw,vec);
OnePhotonView(ψ) == vec
bw = WaveguideBasis(1,3,times);
ψ = onephoton(bw,2,x->dt);
OnePhotonView(ψ,2) == ones(length(times))
source
WaveguideQED.onephotonMethod
onephoton(b::WaveguideBasis{T,Nw},i,ξvec;norm=false) where {T,Nw}
  • Since b contains a Nw waveguides, the index of the waveguide needed.
  • i is the index of the waveguide at which the onephoton wavepacket is created in
  • ξ should be AbstractArray with length(ξ) == b.nsteps
  • norm::Bool=false: If true normalize the resulting wavepacket.
source
WaveguideQED.onephotonMethod
onephoton(b::WaveguideBasis{T,Nw},i::Int,ξ::Function,args...; norm=false) where {T,Nw}
  • Since b contains Nw waveguides, the index of the waveguide needed.
  • i is the index of the waveguide at which the onephoton wavepacket is created in
  • ξ should be broadcastable as ξ.(times,args...), where times = 0:b.dt:(b.nsteps-1)*b.dt (the times used to define the waveguidebasis)
  • args...: additional arguments to be passed to ξ if it is a function.
  • norm::Bool=false: If true normalize the resulting wavepacket.
source
WaveguideQED.onephotonMethod
onephoton(b::WaveguideBasis{T,1},ξ::AbstractArray;norm=false) where {T}
  • Since b only contains a single waveguide, the index of the waveguide is not needed.
  • ξ should be AbstractArray with length(ξ) == b.nsteps
  • norm::Bool=false: If true normalize the resulting wavepacket.
source
WaveguideQED.onephotonMethod
onephoton(b::WaveguideBasis{T,1},ξ::Function,args...,norm=false) where {T}
  • Since b only contains a single waveguide, the index of the waveguide is not needed.
  • ξ should be broadcastable as ξ.(times,args...), where times = 0:b.dt:(b.nsteps-1)*b.dt (the times used to define the waveguidebasis)
  • args...: additional arguments to be passed to ξ if it is a function.
  • norm::Bool=false: If true normalize the resulting wavepacket.
source
WaveguideQED.plot_twophoton!Method
plot_twophoton!(ax,twophotonstate::TwophotonView,times;kwargs...)
plot_twophoton!(ax,twophotonstate::TwoWaveguideView,times;kwargs...)
plot_twophoton!(ax,state::Ket,times;kwargs...)
plot_twophoton!(ax,twophotonstate,times;kwargs...)

Plots the twophoton state in the given ax.

Arguments

  • ax of type PyObject <AxesSubplot: > from PyPlot
  • State to be plotted twophotonstate or state. If state is a Ket TwoPhotonView is called to extract twophotonstate. Otherwise twophotonstate should be AbstractArray of dimensions (length(times),length(times)).

#Return ax.contour object with the plotted state.

source
WaveguideQED.twophotonFunction

Create a twophoton wavepacket of the form $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{j}}^\dagger(t') |\emptyset \rangle$ where $i$ and $j$ are the indeces of the waveguide and return it as a Ket.

See also WaveguideBasis and TwoPhotonView for how to view the state.

#Examples Creating twophoton state with only one waveguide (using a function):

using LinearAlgebra; #hide
times = 1:1:10;
bw = WaveguideBasis(2,times);
ψ = twophoton(bw,(t1,t2)->1,norm=false);
ψview = TwoPhotonView(ψ);
ψview == ones((length(times),length(times)))

ψ = twophoton(bw,(t1,t2,arg1)->arg1,123,norm=false);
ψview = TwoPhotonView(ψ);
ψview == 123*ones((length(times),length(times)))

Creating twophoton state with only one waveguide (using a matrix):

ψ = twophoton(bw,ones((length(times),length(times))),norm=false);
ψview = TwoPhotonView(ψ);
ψview == ones((length(times),length(times)))

Creating twophoton state in waveguide 2 with multiple waveguides

bw = WaveguideBasis(2,3,times)
ψ = twophoton(bw,2,(t1,t2)->1,norm=false);

Creating twophoton state across waveguide 1 and 2

bw = WaveguideBasis(2,3,times)
ψ = twophoton(bw,1,2,(t1,t2)->1,norm=false);
source
WaveguideQED.twophotonMethod
twophoton(b::WaveguideBasis{T,Nw},idx::I,ξ::Function,times,args...;norm=false) where {T,Nw,I<:Union{Vector{Int64},Tuple{Vararg{Int64}}}} = twophoton(b,idx[1],idx[2],ξ,times,args...,norm=norm)

Arguments

  • b::WaveguideBasis{T,Nw} contains multiple waveguides and index i is need.
  • I[1] denotes the index $i$ of the waveguide in the twophoton state: $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{j}}^\dagger(t') |\emptyset \rangle$.
  • I[2] denotes the index $j$ of the waveguide in the twophoton state: $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{j}}^\dagger(t') |\emptyset \rangle$.
  • ξ given as a function should follow ξ(times[l],times[m],args...).
  • times: A vector or range of length(times) = b.nsteps.
  • args...: additional arguments to be passed to ξ.
  • norm::Bool=false: Toggles whether to normalize the resulting wavepacket.
source
WaveguideQED.twophotonMethod
twophoton(b::WaveguideBasis{T,Nw},idx::I,ξ::Matrix;norm=false) where {T,Nw,I<:Union{Vector{Int64},Tuple{Vararg{Int64}}}}

Arguments

  • b::WaveguideBasis{T,Nw} contains multiple waveguides and index i is need.
  • I[1] denotes the index $i$ of the waveguide in the twophoton state: $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{j}}^\dagger(t') |\emptyset \rangle$.
  • I[2] denotes the index $j$ of the waveguide in the twophoton state: $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{j}}^\dagger(t') |\emptyset \rangle$.
  • ξ given as a matrix of dimension (b.nsteps, b.nsteps).
  • norm::Bool=false: Toggles whether to normalize the resulting wavepacket.
source
WaveguideQED.twophotonMethod
twophoton(b::WaveguideBasis{T,Nw},i::Int,ξ::Function,args...;norm=false) where {T,Nw}

Arguments

  • b::WaveguideBasis{T,Nw} contains multiple waveguides and index i is need.
  • i denotes the index of the waveguide in the twophoton state: $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{i}}^\dagger(t') |\emptyset \rangle$.
  • ξ given as a function should follow ξ(t1,t2,args...) where t1 and t2 are the input times.
  • args...: additional arguments to be passed to ξ.
  • norm::Bool=false: Toggles whether to normalize the resulting wavepacket.
source
WaveguideQED.twophotonMethod
twophoton(b::WaveguideBasis{T,Nw},i::Int,j::Int,ξ::Function,args...;norm=false) where {T,Nw}

Arguments

  • b::WaveguideBasis{T,Nw} contains multiple waveguides and index i is need.
  • i denotes the index of the waveguide in the twophoton state: $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{j}}^\dagger(t') |\emptyset \rangle$.
  • j denotes the index of the waveguide in the twophoton state: $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{j}}^\dagger(t') |\emptyset \rangle$.
  • ξ given as a function should follow ξ(t1,t2,args...) where t1 and t2 are the input times.
  • args...: additional arguments to be passed to ξ.
  • norm::Bool=false: Toggles whether to normalize the resulting wavepacket.
source
WaveguideQED.twophotonMethod
twophoton(b::WaveguideBasis{T,Nw},i::Int,j::Int,ξ::Matrix;norm=false) where {T,Nw}

Arguments

  • b::WaveguideBasis{T,Nw} contains multiple waveguides and index i is need.
  • i denotes the index of the waveguide in the twophoton state: $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{j}}^\dagger(t') |\emptyset \rangle$.
  • j denotes the index of the waveguide in the twophoton state: $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{j}}^\dagger(t') |\emptyset \rangle$.
  • ξ given as a matrix of dimension (b.nsteps, b.nsteps).
  • norm::Bool=false: Toggles whether to normalize the resulting wavepacket.
source
WaveguideQED.twophotonMethod
twophoton(b::WaveguideBasis{T,Nw},i::Int,ξ::Matrix;norm=false) where {T,Nw}

Arguments

  • b::WaveguideBasis{T,Nw} contains multiple waveguides and index i is need.
  • i denotes the index of the waveguide in the twophoton state: $\int_{t_0}^{t_{end}} dt' \int_{t_0}^{t_{end}} dt \xi(t,t') w_{\mathrm{i}}^\dagger(t)w_{\mathrm{i}}^\dagger(t') |\emptyset \rangle$.
  • ξ given as a matrix of dimension (b.nsteps, b.nsteps).
  • norm::Bool=false: Toggles whether to normalize the resulting wavepacket.
source
WaveguideQED.twophotonMethod
twophoton(b::WaveguideBasis{T,1},ξ::Function,args...,norm=false) where {T}

Arguments

  • b::WaveguideBasis{T,Nw} contains only one waveguide and no index needed.
  • ξ given as a function should follow ξ(t1,t2,args...), where t1 and t2 are the input times.
  • args...: additional arguments to be passed to ξ if it is a function.
  • norm::Bool=false: Toggles whether to normalize the resulting wavepacket.
source
WaveguideQED.twophotonMethod
twophoton(b::WaveguideBasis{T,1},ξ::Matrix;norm=false) where {T}

Arguments

  • b::WaveguideBasis{T,Nw} contains only one waveguide and no index needed.
  • ξ given as a matrix of dimension (b.nsteps, b.nsteps).
  • norm::Bool=false: Toggles whether to normalize the resulting wavepacket.
source
WaveguideQED.view_waveguideMethod
view_waveguide(ψ::ket)
view_waveguide(ψ::ket,index)

View the Waveguide state given a state ψ containing a WaveguideBasis by returning view(reshape(ψ.data,Tuple(ψ.basis.shape)),index...). If no index is provided the ground state is returned. The index provided should be of the form [:,i,j] where (:) is at the location of the WaveguideBasis and i and j are indeces of other basises. See example:

using QuantumOptics;
times=0:1:10;
bw = WaveguideBasis(1,times);
bc = FockBasis(2);
ψ_waveguide = Ket(bw,ones(length(times)));
ψ_total = ψ_waveguide ⊗ fockstate(bc,0) ⊗ fockstate(bc,0);
ψ_view = view_waveguide(ψ_total);
ψ_view_index = view_waveguide(ψ_total,[:,1,1]);
ψ_view==ψ_view_index
ψ_total = ψ_waveguide ⊗ fockstate(bc,2) ⊗ fockstate(bc,1);
view_waveguide(ψ_total,[:,3,2]) == ones(length(times))
true
source
WaveguideQED.waveguide_evolutionMethod
waveguide_evolution(tspan, psi0, H; fout)

Integrate time-dependent Schroedinger equation to evolve states or compute propagators.

Arguments

  • tspan: Vector specifying the points of time for which output should be displayed.
  • psi0: Initial state vector can only be a ket.
  • H: Operator containing a WaveguideOperator either through a LazySum or LazyTensor.
  • fout=nothing: If given, this function fout(t, psi) is called every time step. Example: fout(t,psi) = expect(A,psi) will return the epectation value of A at everytimestep. If fout =1 the state psi is returned for all timesteps in a vector. ATTENTION: The state psi is neither normalized nor permanent! It is still in use by the ode solver and therefore must not be changed.

Returns

  • if fout=nothing the output of the solver will be the state ψ at the last timestep.
  • if fout is given a tuple with the state ψ at the last timestep and the output of fout is given. If fout returns a tuple the tuple will be flattened.
  • if fout = 1 ψ at all timesteps is returned.

Examples

  • fout(t,psi) = (expect(A,psi),expect(B,psi)) will result in a tuple (ψ, ⟨A(t)⟩,⟨B(t)⟩), where ⟨A(t)⟩ is a vector with the expectation value of A as a function of time.
source
LinearAlgebra.mul!Function
mul!(result::Ket{B1}, a::LazyTensor{B1,B2,F,I,T}, b::Ket{B2}, alpha, beta)

In-place multiplication of operators/state vectors. Updates result as result = alpha*a*b + beta*result. a is a LazyTensor that contains a WaveguideOperator

source
mul!(result::Ket{B1}, a::CavityWaveguideEmission, b::Ket{B2}, alpha, beta) where {B1<:Basis,B2<:Basis}
mul!(result::Ket{B1}, a::CavityWaveguideAbsorption, b::Ket{B2}, alpha, beta) where {B1<:Basis,B2<:Basis}

Fast in-place multiplication of operators/state vectors. Updates result as result = alpha*a*b + beta*result. a is a CavityWaveguideOperator.

source
QuantumOpticsBase.createFunction
create(basis::WaveguideBasis{Np,1}) where {Np}
create(basis::WaveguideBasis{Np,Nw},i::Int) where {Np,Nw}

Creation operator $w^\dagger$ for WaveguideBasis $w_i(t_k)^\dagger | \emptyset \rangle = | 1_k \emptyset \rangle_i$ with cutoff (maximum number of photons Np) where i is the index of the waveguide. $t_k$ is determined by the timeindex property of the operator which can be changed by set_waveguidetimeindex!(a::WaveguideOperator,k::Int)

Arguments

  • basis of type WaveguideBasis, defines cutoff photon number Np and number of waveguides Nw
  • i determines which waveguide the operator acts on and should be i ≤ Nw. If Nw=1 then i=1 is assumed (there is only one waveguide).

Returns

WaveguideCreate

source
QuantumOpticsBase.destroyFunction
destroy(basis::WaveguideBasis{Np,1}) where {Np}
destroy(basis::WaveguideBasis{Np,Nw},i::Int) where {Np,Nw}

Annihilation operator $w$ for WaveguideBasis $w_i(t_k) | 1_j \emptyset \rangle_i = \delta_{k,j} | \emptyset \rangle$ with cutoff (maximum number of photons Np) where i is the index of the waveguide. $t_k$ is determined by the timeindex property of the operator which can be changed by set_waveguidetimeindex!(a::WaveguideOperator,k::Int)

Arguments

  • basis of type WaveguideBasis, defines cutoff photon number Np and number of waveguides Nw
  • i determines which waveguide the operator acts on and should be i ≤ Nw. If Nw=1 then i=1 is assumed (there is only one waveguide).

Returns

WaveguideDestroy

source
QuantumInterface.tensorFunction
tensor(a::AbstractOperator,b::CavityWaveguideAbsorption)
tensor(a::CavityWaveguideAbsorption,b::AbstractOperator)
tensor(a::AbstractOperator,b::CavityWaveguideEmission)
tensor(a::CavityWaveguideEmission,b::AbstractOperator)

Methods for tensorproducts between QuantumOptics.jl operators and CavityWaveguideOperator.

source