Processing.jl

Luna.Processing.AutoWindowType
AutoWindow(width, λmin, λmax, ω0fun; relative=false, ndims=1)

Window function generator which automatically tracks the central frequency in the spectral region given by λmin and λmax and applies a window of a specific width around the peak. The central frequency is found using the function ω0fun(ω, Iω::AbstractVector), where ω and are already cropped to within the wavelength limits given. If relative is true, width is relative bandwidth instead of the wavelength width. ndims determines how many dimensions of the array to sum over. For a field array with size (Nω, N1, N2, ...), the first dimension is always assumed to be frequency. ndim=1 means each field to be analysed is 1-dimensional, so the window iterates over all of (N1, N2, ...). ndim=2 means each field to be analysed is 2-dimensional, (Nω, N1) in size, and will be summed over its second dimension before finding the central frequency. The window iterates over all other dimensions, (N2, ...).

A AutoWindow automatically stores the limits of the windows it applies in the field lims.

source
Luna.Processing.CommonType
Common(val)

Wrapper type to tell scanproc that val is the same for each simulation being processed, and so only needs to be returned once rather than for each simulation in the scan.

source
Luna.Processing.VarLengthType
VarLength(val)

Wrapper type to tell scanproc that the shape of val is different for each simulation being processed. Return values wrapped in VarLength will be placed in an array of arrays.

Note

While the shape of val can be different between simulations, the type must be the same, including the dimensionality and element type of arrays.

source
Luna.Processing.CentroidWindowMethod
CentroidWindow(width, λmin, λmax; relative=false, ndims=1, power=1)

An AutoWindow which uses the centroid (centre of mass or first moment) of the spectral energy density as the central frequency. Before calculating the centroid, the SED is raised to the power given.

source
Luna.Processing._specres_kernel!Method

Convolution kernel for each output point. We simply loop over all outer indices and output points. The inner loop adds up the contributions from the specified window around the target point. Note that this works without scaling also for wavelength ranges because the integral is still over a frequency grid (with appropriate frequency dependent integration bounds).

source
Luna.Processing.arrivaltimeMethod
arrivaltime(grid, Eω; bandpass=nothing, method=:moment, oversampling=1)

Extract the arrival time of the pulse in the wavelength limits λlims.

Arguments

  • bandpass : method to bandpass the field if required. See window_maybe
  • method::Symbol : :moment to use 1st moment to extract arrival time, :peak to use the time of peak power
  • oversampling::Int : If >1, oversample the time-domain field before extracting delay
  • sumdims : Single Int or Tuple of Ints. The time-domain power will be summed over these dimensions (e.g. modes) before extracting the arrival time.
source
Luna.Processing.beamMethod
beam(grid, Eωm, modes, x, y; z=0, components=:xy)
beam(output, x, y, zslice; bandpass=nothing)

Calculate the beam profile of the multi-mode field Eωm on the grid given by spatial coordinates x and y. If output is given, create the modes from that and take the field nearest propagation slice zslice.

source
Luna.Processing.coherenceMethod
coherence(Eω; ndim=1)

Calculate the first-order coherence function g₁₂ of the set of fields . The ensemble average is taken over the last ndim dimensions of , other dimensions are preserved.

See J. M. Dudley and S. Coen, Optics Letters 27, 1180 (2002).

source
Luna.Processing.envelopeMethod
envelope(grid, Eω)

Get the envelope electric field including the carrier wave from the frequency-domain field sampled on grid.

source
Luna.Processing.fwhm_fMethod
fwhm_f(grid, Eω::Vector; bandpass=nothing, oversampling=1, sumdims=nothing, minmax=:min)

Extract the frequency FWHM. If bandpass is given, bandpass the field according to window_maybe. If sumdims is given, the energy density is summed over these dimensions (e.g. modes) before extracting the FWHM. minmax determines determines whether the FWHM is taken at the narrowest (:min) or the widest (:max) point.

source
Luna.Processing.fwhm_tMethod
fwhm_t(grid::AbstractGrid, Eω; bandpass=nothing, oversampling=1, sumdims=nothing, minmax=:min)

Extract the temporal FWHM. If bandpass is given, bandpass the fieldaccording to window_maybe. If oversampling > 1, the time-domain field is oversampled before extracting the FWHM. If sumdims is given, the time-domain power is summed over these dimensions (e.g. modes) before extracting the FWHM. minmax determines determines whether the FWHM is taken at the narrowest (:min) or the widest (:max) point.

source
Luna.Processing.getEtMethod
getEt(grid, Eω; trange=nothing, oversampling=4, bandpass=nothing, FTL=false)

Get the envelope time-domain electric field (including the carrier wave) from the frequency- domain field . The field can be cropped in time using trange, it is oversampled by a factor of oversampling (default 4) and can be bandpassed with bandpass (see window_maybe). If FTL is true, return the Fourier-transform limited pulse, i.e. remove any spectral phase.

If zslice is given, returs only the slices of closest to the given distances. zslice can be a single number or an array.

source
Luna.Processing.getEtMethod
getEt(output[, zslice]; kwargs...)

Get the envelope time-domain electric field (including the carrier wave) from the output. If zslice is given, returs only the slices of closest to the given distances. zslice can be a single number or an array.

source
Luna.Processing.getEtxyMethod
getEtxy(output, xs, z; kwargs...)
getEtxy(Etm, modes, xs, z; components=:xy)

Calculate the time-dependent electric field at transverse position xs and longitudinal position z from either the modal time-dependent field Etm or the given output.

xs should be a 2-Tuple of coordinates–either (r, θ) for polar coordinates or (x, y) in Cartesian coordinates, depending on the coordinate system of the modes–or a 2-Tuple of vectors containing the coordinates. If vectors are given, the output contains values of Etxy at all combinations of the coordinates.

Additional keyword arguments to getEtxy(output, ...) are passed through to Processing.getEt

source
Luna.Processing.getEωMethod
getEω(output[, zslice])

Get frequency-domain modal field from output with correct normalisation (i.e. abs2.(Eω)` gives angular-frequency spectral energy density in J/(rad/s)).

source
Luna.Processing.getIωMethod
getIω(ω, Eω, specaxis; specrange=nothing, resolution=nothing)

Get spectral energy density and x-axis given a frequency array ω and frequency-domain field , assumed to be correctly normalised (see getEω). specaxis determines the x-axis:

  • :f -> x-axis is frequency in Hz and Iω is in J/Hz
  • :ω -> x-axis is angular frequency in rad/s and Iω is in J/(rad/s)
  • :λ -> x-axis is wavelength in m and Iω is in J/m

Keyword arguments

  • specrange::Tuple can be set to a pair of limits on the spectral range (in specaxis units).
  • resolution::Real is set, smooth the spectral energy density as defined by specres.

Note that if resolution and specaxis=:λ is set it is highly recommended to also set specrange.

source
Luna.Processing.getIωMethod
getIω(output, specaxis[, zslice]; kwargs...)

Calculate the correctly normalised frequency-domain field and convert it to spectral energy density on x-axis specaxis (:f, , or ). If zslice is given, returs only the slices of closest to the given distances. zslice can be a single number or an array. specaxis determines the x-axis:

  • :f -> x-axis is frequency in Hz and Iω is in J/Hz
  • :ω -> x-axis is angular frequency in rad/s and Iω is in J/(rad/s)
  • :λ -> x-axis is wavelength in m and Iω is in J/m

Keyword arguments

  • specrange::Tuple can be set to a pair of limits on the spectral range (in specaxis units).
  • resolution::Real is set, smooth the spectral energy density as defined by specres.

Note that resolution is set and specaxis=:λ it is highly recommended to also set specrange.

source
Luna.Processing.ionisation_fractionMethod
ionisation_fraction(output, xs; ratefun, oversampling=1)
ionisation_fraction(output; ratefun, oversampling=1, maxevals=1000)

Calculate the ionisation fraction at transverse coordinates xs using the ionisation-rate function ratefun. If xs is not given, calculate the average ionisation fraction across the waveguide core. In this case, maxevals determines the maximum number of function evaluations for the integral.

Warning

Calculating the average ionisation fraction is much slower than calculating it at a single point

source
Luna.Processing.makemodesMethod
makemodes(output)

Create the modes used in a simulation using MarcatiliModes. If output was created by Interface.prop_capillary_args and hence has a field prop_capillary_args, this is used to match the gas fill from the simulation. Otherwise, the modes are created without gas fill.

source
Luna.Processing.nearest_zMethod
nearest_z(output, z)

Return the index of saved z-position(s) closest to the position(s) z. Output is always an array, even if z is a number. If z is negative, its absolute value is taken as the fraction of the total propagation distance.

source
Luna.Processing.peakpowerMethod
peakpower(grid, Eω; bandpass=nothing, oversampling=1, sumdims=nothing)
peakpower(output; bandpass=nothing, oversampling=1, sumdims=nothing)

Extract the peak power. If bandpass is given, bandpass the field according to window_maybe. If sumdims is not nothing, sum the time-dependent power over these dimensions (e.g. modes) before taking the maximum.

source
Luna.Processing.scanprocMethod
scanproc(f, scanfiles)
scanproc(f, directory)
scanproc(f, directory, pattern)
scanproc(f)

Iterate over the scan output files, apply the processing function f(o::AbstractOutput), and collect the results in arrays.

The files can be given as:

  • a Vector of AbstractStrings containing file paths
  • a directory to search for files according to the naming pattern of Output.ScanHDF5Output
  • a directory and a glob pattern

If nothing is specified, scanproc uses the current working directory.

f can return a single value, an array, or a tuple/array of arrays/numbers. Arrays returned by f must either be of the same size for each processed file, or wrapped in a VarLength. Values returned by f which are guaranteed to be identical for each processed file can be wrapped in a Common, and scanproc only returns these once.

Example

Et, Eω = scanproc("path/to/scandir") do output
    t, Et = getEt(output)
    ω, Eω = getEω(output)
    energyout = energyout = Processing.VarLength(output["stats"]["energy"])
    Common(t), Et, Common(ω), Eω, energyout
end
source
Luna.Processing.scanprocMethod
scanproc(f, outputs; shape=nothing)

Iterate over the scan outputs, apply the processing function f(o::AbstractOutput), and collect the results in arrays.

If the outputs are MemoryOutputs which do not contain the scan metadata, the shape of the scan must be given explicitly (e.g. via size(scan)).

f can return a single value, an array, or a tuple/array of arrays/numbers. Arrays returned by f must either be of the same size for each processed output, or wrapped in a VarLength. Values returned by f which are guaranteed to be identical for each processed output can be wrapped in a Common, and scanproc only returns these once.

source
Luna.Processing.specresMethod
specres(ω, Iω, specaxis, resolution, specrange; window=nothing, nsamples=10)

Smooth the spectral energy density Iω(ω) to account for the given resolution on the defined specaxis and specrange. The window function to use defaults to a Gaussian function with FWHM of resolution, and by default we sample nsamples=10 times within each resolution.

Note that you should prefer the resolution keyword of getIω instead of calling this function directly.

The input ω and should be as returned by getIω with specaxis = :ω.

Returns the new specaxis grid and smoothed spectrum.

source
Luna.Processing.time_bandwidthMethod
time_bandwidth(grid, Eω; bandpass=nothing, oversampling=1)

Extract the time-bandwidth product, after bandpassing if required. The TBP is defined here as ΔfΔt where Δx is the FWHM of x. (In this definition, the TBP of a perfect Gaussian pulse is ≈0.44). If oversampling > 1, the time-domain field is oversampled before extracting the FWHM.

source
Luna.Processing.window_maybeMethod
window_maybe(ω, Eω, win)

Apply a frequency window to the field if required. Possible values for win:

  • nothing : no window is applied
  • 4-Tuple of Numbers : the 4 parameters for a Maths.planck_taper in wavelength
  • 3-Tuple of Numbers : minimum, maximum wavelength, and smoothing in radial frequency
  • 2-Tuple of Numbers : minimum and maximum wavelength with automatically chosen smoothing
  • Vector{<:Real} : a pre-defined window function (shape must match ω)
  • PeakWindow : automatically track the peak in a given range and apply the window around it
  • window(ω, Eω) : an arbitrary user-supplied window function
source
Luna.Processing.ωwindow_λMethod
ωwindow_λ(ω, λlims; winwidth=:auto)

Create a ω-axis filtering window to filter in λlims. winwidth, if a Number, sets the smoothing width of the window in rad/s.

source