Technical reference

Public interface

Peacock.BrillouinZoneCoordinateType
BrillouinZoneCoordinate(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.

source
Peacock.GeometryType
Geometry(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.

source
Peacock.GeometryMethod
Geometry(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.

source
Peacock.GeometryMethod
Geometry(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.

source
Peacock.HilbertSpaceType
HilbertSpace(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.

source
Peacock.HilbertSpaceMethod
HilbertSpace(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.

source
Peacock.ModeType
Mode(k0, frequency, data, weighting, basis, label)

Eigenmode of a photonic crystal expressed on a plane-wave basis with a weighted inner product.

source
Peacock.SolverType
Solver(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.

source
Peacock.SolverMethod
Solver(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.

source
Peacock.SolverMethod
Solver(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).

source
Peacock.SolverMethod
Solver(geometry::Geometry, basis::PlaneWaveBasis)

Approximate the geometry using the given basis of plane waves.

source
Peacock.get_fieldMethod
get_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.

source
Peacock.plot_band_diagramMethod
plot_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 points
  • labels=[]: overwrite the labels for each k in ks
  • bands=(:): indices of the bands to plot
  • frequency_scale=1: rescales the frequencies before plotting
  • color="k": color of the bands
  • markersize=nothing: overwrite the size of each point
  • show_vlines=true: plot vertical lines at each k in ks
source
Peacock.plot_wilson_loop_windingMethod
plot_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 same k in different Brillouin zones
  • labels=[]: overwrite the labels for each k in ks
  • markersize=nothing: overwrite the size of each point
source
Peacock.solveFunction
solve(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.

source
Peacock.solveMethod
solve(solver::Solver, k::AbstractVector{<:Real}, polarisation::Polarisation; bands=:)

Calculate the eigenmodes of a photonic crystal at position k in reciprocal space.

source
PyPlot.plotMethod
plot(geometry::Geometry)

Plot the permittivity and permeability of the geometry in real space.

source
PyPlot.plotMethod
plot(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.

source
PyPlot.plotMethod
plot(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.

source

Private interface

The `Peacock.jl` namespace

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.DiagonalMatrixType
DiagonalMatrix(diag::AbstractVector{ComplexF64})

A sparse diagonal matrix that can be used in left division (D \ X)

source
Peacock.PlaneWaveBasisType
PlaneWaveBasis(b1::Vector{Float64}, b2::Vector{Float64},
        ps::Vector{Int}, qs::Vector{Int})

A 2D basis of plane waves, ks = ps*b1 + qs*b2.

source
Peacock.PlaneWaveBasisMethod
Solver(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.

source
Peacock.PlaneWaveBasisMethod
PlaneWaveBasis(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).

source
Peacock.as_to_bsMethod
as_to_bs(a1, a2)

Calculate reciprocal lattice vectors, b1 and b2, from the real space lattice vectors, a1 and a2.

source
Peacock.bs_to_asMethod
bs_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.

source
Peacock.convmatMethod
convmat(mat::AbstractMatrix, basis::PlaneWaveBasis)

Generate convolution matrices, see Raymond Rumpf's CEM Lecture #18, "Maxwell's Equations in Fourier Space", for further reading.

source
Peacock.get_kMethod
get_k(coord::BrillouinZoneCoordinate, basis::PlaneWaveBasis)

Return the k-space coordinate of the BrillouinZoneCoordinate in a particular PlaneWaveBasis, ie k = p*b1 + q*b2.

source
Peacock.normaliseMethod
normalise(data; weighting=I)

Normalisation of vectors with a weighted inner product.

source
Peacock.orthonormaliseMethod
orthonormalise(data; weighting=I)

Gram-Schmidt orthonormalisation of vectors with a weighted inner product.

source
Peacock.overlapsMethod
overlaps(a::HilbertSpace, b::HilbertSpace)

Calculate the overlaps between the basis vectors of each Hilbert space.

source
Peacock.plot_fieldMethod
plot_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.

source
Peacock.plot_fieldMethod
plot_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.

source
Peacock.plot_fieldMethod
plot_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.

source
Peacock.sample_pathMethod
sample_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].

source
Peacock.shift_k0Method
shift_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.

source
Peacock.unitary_approxMethod
unitary_approx(M::AbstractArray)

Calculate the best unitary approxmation of M using singular value decomposition.

source
Peacock.wilson_eigenMethod
wilson_eigvals(spaces::AbstractArray{HilbertSpace,1}; closed=true)

Return the eigenvalues and eigenvectors of wilson_matrix, sorted by the phase angle of the eigenvalues.

source
Peacock.wilson_gaugeMethod
wilson_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.

source
Peacock.wilson_matrixMethod
wilson_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.

source