Technical reference
Public interface
Peacock.BrillouinZoneCoordinate
— TypeBrillouinZoneCoordinate(p::Float64, q::Float64, label::String="")
A labelled coordinate in the Brillouin zone.
The arguments p
and q
are the coefficients of reciprocal lattice vectors b1
and b2
. The k-space coordinate, k = p * b1 + q * b2
, is generated by get_k(coord::BrillouinZoneCoordinate, basis::PlaneWaveBasis)
. For example, BrillouinZoneCoordinate(0.5,0)
is on the edge of the first Brillouin zone.
Peacock.Geometry
— TypeGeometry(a1::Array{Real,1}, a2::Array{Real,1}, ep::Array{ComplexF64,2}, mu::Array{ComplexF64,2})
2D geometry defined in real space, with lattice vectors a1
and a2
, and relative permeability and permittivity ep
and mu
, respectively.
Peacock.Geometry
— MethodGeometry(epf::Function, muf::Function, a1::Array{<:Real,1}, a2::Array{<:Real,1}, d1::Real, d2::Real)
Generate geometry with permittivity epf(x,y)
and permeability muf(x,y)
.
The real space lattice vectors, a1
and a2
, define the unit cell. The grid resolution along each lattice vector is d1
and d2
, respectively.
Peacock.Geometry
— MethodGeometry(epf::Function, muf::Function, a1::Array{<:Real,1}, a2::Array{<:Real,1}, d1::Real, d2::Real)
Generate geometry with permittivity epf(x,y)
and permeability muf(x,y)
.
The real space lattice vectors are assumed to have unit length and are at angles a1_deg
and a2_deg
counter-clockwise from the x-axis. The grid resolution along each lattice vector is d1
and d2
, respectively.
Peacock.HilbertSpace
— TypeHilbertSpace(k0, data, weighting, basis)
Hilbert space spanned by the eigenvectors in each column of data
.
The eigenvectors will be orthonormalised using Gram-Schmidt orthonormalisation, see orthonormalise
.
Peacock.HilbertSpace
— MethodHilbertSpace(modes::Array{Mode,1})
Returns the Hilbert space spanned by the modes
.
The data
of the Hilbert space is guaranteed to be orthonormal under the weighting of the modes
.
Peacock.Mode
— TypeMode(k0, frequency, data, weighting, basis, label)
Eigenmode of a photonic crystal expressed on a plane-wave basis
with a weighted inner product.
Peacock.Solver
— TypeSolver(basis::PlaneWaveBasis, epc::Matrix{ComplexF64}, muc::Matrix{ComplexF64})
A plane-wave expansion method solver, where epc
and muc
are the convolution matrices of the permittivity and permeability, respectively, for the given basis
of plane waves.
Peacock.Solver
— MethodSolver(geometry::Geometry, cutoff_b1::Int, cutoff_b2::Int)
Approximate the geometry using a basis of plane waves truncated in a rhombus.
The rhombus has lengths cutoff_b1
and cutoff_b2
in the b1
and b2
directions, respectively.
Peacock.Solver
— MethodSolver(geometry::Geometry, cutoff::Int)
Approximate the geometry using a basis of plane waves truncated in a circle.
The circle has a diameter of cutoff
Brillouin zones. Increasing the cutoff
will increase the number of plane waves leading to a more accurate solution. It is assumed that norm(b1) == norm(b2)
.
Peacock.Solver
— MethodSolver(geometry::Geometry, basis::PlaneWaveBasis)
Approximate the geometry using the given basis
of plane waves.
Peacock.get_field
— Methodget_field(data::AbstractVector{<:Complex}, basis::PlaneWaveBasis;
k0=[0,0], t1s=-0.5:0.01:0.5, t2s=-0.5:0.01:0.5)
Convert the data
from a PlaneWaveBasis
to a real space grid.
Peacock.plot_band_diagram
— Methodplot_band_diagram(my_solve::Function, ks, <keyword arguments>)
Plot the bands generated by my_solve(k)
along ks
.
Keyword arguments
dk=nothing
: maximum distance between pointslabels=[]
: overwrite the labels for eachk
inks
bands=(:)
: indices of the bands to plotfrequency_scale=1
: rescales the frequencies before plottingcolor="k"
: color of the bandsmarkersize=nothing
: overwrite the size of each pointshow_vlines=true
: plot vertical lines at eachk
inks
Peacock.plot_band_diagram
— Methodplot_band_diagram(solver::Solver, ks, polarisation::Polarisation, <keyword arguments>)
Plot the bands generated by solve(solver, k, polarisation)
along ks
.
Takes the same keyword arguments as plot_band_diagram(my_solve::Function, ks)
, but here ks
can also include BrillouinZoneCoordinate
s.
Peacock.plot_wilson_loop_winding
— Methodplot_wilson_loop_winding(solver::Solver, ks, polarisation, bands::AbstractVector{<:Int}, <keyword arguments>)
Perform a series of Wilson loops along ks
, and plot the spectra on a band diagram.
Keyword arguments
dk_outer=nothing
: maximum distance between each loop (resolution of the scan)dk_inner=nothing
: maximum distance between points along a loop (resolution of the loop)delta_brillouin_zone=BrillouinZoneCoordinate(0,1)
: each Wilson loop starts at and finishes in at the samek
in different Brillouin zoneslabels=[]
: overwrite the labels for eachk
inks
markersize=nothing
: overwrite the size of each point
Peacock.solve
— Functionsolve(solver::Solver, x::BrillouinZoneCoordinate, polarisation::Polarisation; bands=:)
Calculate the eigenmodes of a photonic crystal at position k=x.p*b1 + x.q*b2
in reciprocal space.
Peacock.solve
— Methodsolve(solver::Solver, k::AbstractVector{<:Real}, polarisation::Polarisation; bands=:)
Calculate the eigenmodes of a photonic crystal at position k
in reciprocal space.
PyPlot.plot
— Methodplot(geometry::Geometry)
Plot the permittivity and permeability of the geometry
in real space.
PyPlot.plot
— Methodplot(mode::Mode, [bloch_phase=true])
Plot the mode in real space.
To plot only the cell-periodic part of the Bloch wave, set bloch_phase=false
.
PyPlot.plot
— Methodplot(solver::Solver)
Plot the representation of the solver
's geometry in real space.
The solver
approximates the geometry using a truncated basis of plane waves, so plot(solver)
lets us judge how accurately the geometry is represented.
Private interface
The following objects are used internally, and aren't exported to the global namespace. However, they may be useful for more advanced users, in which case you should use Peacock.[name]
to access these.
Peacock.DiagonalMatrix
— TypeDiagonalMatrix(diag::AbstractVector{ComplexF64})
A sparse diagonal matrix that can be used in left division (D \ X)
Peacock.PlaneWaveBasis
— TypePlaneWaveBasis(b1::Vector{Float64}, b2::Vector{Float64},
ps::Vector{Int}, qs::Vector{Int})
A 2D basis of plane waves, ks = ps*b1 + qs*b2
.
Peacock.PlaneWaveBasis
— MethodSolver(geometry::Geometry, cutoff_b1::Int, cutoff_b2::Int)
Approximate the geometry using a basis of plane waves truncated in a rhombus.
The rhombus has lengths cutoff_b1
and cutoff_b2
in the b1
and b2
directions, respectively.
Peacock.PlaneWaveBasis
— MethodPlaneWaveBasis(geometry::Geometry, cutoff::Int)
Approximate a basis of plane waves truncated in a circle.
The circle has a diameter of cutoff
Brillouin zones. Increasing the cutoff
will increase the number of plane waves leading to a more accurate solution. It is assumed that norm(b1) == norm(b2)
.
Peacock.Polarisation
— TypeTransverse electric (TE) or transverse magnetic (TM) polarisation.
Peacock.as_to_bs
— Methodas_to_bs(a1, a2)
Calculate reciprocal lattice vectors, b1
and b2
, from the real space lattice vectors, a1
and a2
.
Peacock.bs_to_as
— Methodbs_to_as(b1, b2)
Calculate real space lattice vectors, a1
and a2
, from reciprocal lattice vectors, b1
and b2
.
This is actually the same as as_to_bs
, but I think having both functions makes the intention of the code more obvious.
Peacock.convmat
— Methodconvmat(mat::AbstractMatrix, basis::PlaneWaveBasis)
Generate convolution matrices, see Raymond Rumpf's CEM Lecture #18, "Maxwell's Equations in Fourier Space", for further reading.
Peacock.get_k
— Methodget_k(coord::BrillouinZoneCoordinate, basis::PlaneWaveBasis)
Return the k-space coordinate of the BrillouinZoneCoordinate
in a particular PlaneWaveBasis
, ie k = p*b1 + q*b2
.
Peacock.normalise
— Methodnormalise(data; weighting=I)
Normalisation of vectors with a weighted inner product.
Peacock.orthonormalise
— Methodorthonormalise(data; weighting=I)
Gram-Schmidt orthonormalisation of vectors with a weighted inner product.
Peacock.overlaps
— Methodoverlaps(a::HilbertSpace, b::HilbertSpace)
Calculate the overlaps between the basis vectors of each Hilbert space.
Peacock.plot_field
— Methodplot_field(data::AbstractVector{<:Complex}, basis::PlaneWaveBasis;
k0=[0,0], cmap="coolwarm", vmin=nothing, vmax=nothing, label=nothing)
Converts the plane-wave amplitudes to real space using get_field
and then plots the field.
Peacock.plot_field
— Methodplot_field(field::Array{<:Real,2}, a1::Array{<:Real,1}, a2::Array{<:Real,1}; cmap="coolwarm", vmin=nothing, vmax=nothing)
Plot the complex-valued field
on a unit cell with lattice vectors a1
and a2
.
Peacock.plot_field
— Methodplot_field(field::Array{<:Real,2}, a1::Array{<:Real,1}, a2::Array{<:Real,1}; cmap="coolwarm", vmin=nothing, vmax=nothing)
Plot the real-valued field
on a unit cell with lattice vectors a1
and a2
.
Peacock.sample_path
— Methodsample_path(ks; [labels=[]], [dk=nothing])
Return a path through ks
where dk
is the maximum distance between points.
If dk==nothing
then there will be approximately 10 points between ks[1] and ks[2].
Peacock.shift_k0
— Methodshift_k0(space::HilbertSpace, dp::Int, dq::Int)
Shift the basis of the Hilbert space by dp*b1 + dq*b2
, where b1
and b2
are reciprocal lattice vectors.
This is required when we need the overlaps of modes that are at the same k-point but in different Brillouin zones.
Peacock.unitary_approx
— Methodunitary_approx(M::AbstractArray)
Calculate the best unitary approxmation of M
using singular value decomposition.
Peacock.unitary_overlaps
— Methodunitary_overlaps(a::HilbertSpace, b::HilbertSpace)
Shortcut for unitary_approx(overlaps(a, b))
Peacock.wilson_eigen
— Methodwilson_eigvals(spaces::AbstractArray{HilbertSpace,1}; closed=true)
Return the eigenvalues and eigenvectors of wilson_matrix
, sorted by the phase angle of the eigenvalues.
Peacock.wilson_eigvals
— Methodwilson_eigvals(spaces::AbstractArray{HilbertSpace,1}; closed=true)
Return the eigenvalues of wilson_matrix
, sorted by phase angle.
Peacock.wilson_gauge
— Methodwilson_gauge(spaces::AbstractArray{HilbertSpace,1}; closed=true)
Return the eigenvalues, eigenvectors, and gauge of the Wilson loop through the Hilbert space, sorted by the phase angle of the eigenvalues.
Peacock.wilson_matrix
— Methodwilson_matrix(spaces::Array{HilbertSpace,1}; closed::Bool=true)
Calculate the Wilson loop matrix through the Hilbert spaces.
The closed
keyword will assert that spaces[1] == spaces[end]
. Otherwise, it will be assumed that the Wilson loop begins and finishes at the same k0
, but in different Brillouin zones.