ElectricFields.BK7Constant
BK7

Borosilicate glass is commonly used in optical lenses

\[n^2(\lambda) = 1 + \frac{1.03961212\lambda^2}{\lambda^2 - 6.00069867\times10^{-3}\;\textrm{μm}^2} + \frac{0.231792344\lambda^2}{\lambda^2 - 2.00179144\times10^{-2}\;\textrm{μm}^2} + \frac{1.01046945}{\lambda^2 - 1.03560653\times10^{2}\;\textrm{μm}^2}.\]

ElectricFields.CalciteConstant
Calcite

Calcite is a uniaxial, birefringent crystal.

Numbers taken from http://www.newlightphotonics.com/v1/calcite-properties.html.

\[\begin{aligned} n_{\mathrm{o}}^2(\lambda) &= 2.69705 + \frac{0.0192064}{\lambda^2-0.01820\;\textrm{μm}^2} - 0.0151624\lambda^2, \\ n_{\mathrm{e}}^2(\lambda) &= 2.18438 + \frac{0.0087309}{\lambda^2-0.01018\;\textrm{μm}^2} - 0.0024411\lambda^2. \end{aligned}\]

ElectricFields.KTPConstant
KTP

KTP is a biaxial crystal.

Numbers taken from https://www.redoptronics.com/KTP-crystal.html.

\[\begin{aligned} n_x^2(\lambda) &= 2.10468 + \frac{0.89342\lambda^2}{\lambda^2-0.04438\;\textrm{μm}^2} - 0.01036\;\textrm{μm}^{-2}\lambda^2, \\ n_y^2(\lambda) &= 2.14559 + \frac{0.87629\lambda^2}{\lambda^2-0.0485\;\textrm{μm}^2} - 0.01173\;\textrm{μm}^{-2}\lambda^2, \\ n_z^2(\lambda) &= 1.9446 + \frac{1.3617\lambda^2}{\lambda^2-0.047\;\textrm{μm}^2} - 0.01491\;\textrm{μm}^{-2}\lambda^2. \\ \end{aligned}\]

ElectricFields.QuartzConstant
Quartz

Quartz is a uniaxial, birefrigent crystal.

Numbers taken from https://www.newlightphotonics.com/v1/quartz-properties.html.

\[\begin{aligned} n_{\mathrm{o}}^2(\lambda) &=2.3573 - 0.01170\;\textrm{μm}^{-2}\lambda^2 + \frac{0.01054\;\textrm{μm}^{2}}{\lambda^2} + \frac{1.3414\times10^{-4}\;\textrm{μm}^{4}}{\lambda^4} - \frac{4.4537\times10^{-7}\;\textrm{μm}^{6}}{\lambda^6} + \frac{5.9236\times10^{-8}\;\textrm{μm}^{8}}{\lambda^8}, \\ n_{\mathrm{e}}^2(\lambda) &=2.3849 - 0.01259\;\textrm{μm}^{-2}\lambda^2 + \frac{0.01079\;\textrm{μm}^{2}}{\lambda^2} + \frac{1.6518\times10^{-4}\;\textrm{μm}^{4}}{\lambda^4} - \frac{1.9474\times10^{-6}\;\textrm{μm}^{6}}{\lambda^6} + \frac{9.3648\times10^{-8}\;\textrm{μm}^{8}}{\lambda^8}. \end{aligned}\]

ElectricFields.SiO₂Constant
SiO₂

Fused silica

\[n^2(\lambda) = 1 + \frac{0.6961663\lambda^2}{\lambda^2-(0.0684043\;\textrm{μm})^2} + \frac{0.4079426\lambda^2}{\lambda^2-(0.1162414\;\textrm{μm})^2} + \frac{0.8974794\lambda^2}{\lambda^2-(9.896161\;\textrm{μm})^2}.\]

ElectricFields.AbstractKnotSetType
AbstractKnotSet{k,ml,mr,T}

Abstract base for B-spline knot sets. T is the eltype of the knot set, k is the order of the piecewise polynomials (order = degree + 1) and ml and mr are the knot multiplicities of the left and right endpoints, respectively.

ElectricFields.ArbitraryKnotSetType
ArbitraryKnotset(k, t[, ml=k, mr=k])

Construct an order-k arbitrary knot set, with the locations of the knots given by non-decreasing vector t.

ElectricFields.ArbitraryKnotSetType
ArbitraryKnotSet{k,ml,mr}(t)

An arbitrary knot set of order k and left and right multiplicities of ml and mr, respectively. The knot set is specified by the non-decreasing vector t; the same knot may appear multiple times, which influences the continuity of the B-splines at that location.

ElectricFields.BSplineType
BSpline(t, x, w, B, S)

Basis structure for the B-splines generated by the knot set t.

ElectricFields.BSplineFieldType
BSplineField(B, C)

Represents an electric field using an expansion over a B-spline basis B, with C being the array (vector for 1d fields, 3-column matrix for 3d fields) expansion coefficients for the vector potential.

ElectricFields.BlackmanType
ElectricFields.Blackman

Blackman is a sum-of-cosines window, with the terms

\[W(x) = 0.42 + 0.5\cos(2 \pi x) + 0.08\cos(4 \pi x).\]

ElectricFields.BlackmanExactType
ElectricFields.BlackmanExact

Blackman exact is a sum-of-cosines window, with the terms

\[W(x) = 0.4265907136715391 + 0.4965606190885641\cos(2 \pi x) + 0.07684866723989682\cos(4 \pi x).\]

ElectricFields.BlackmanHarrisType
ElectricFields.BlackmanHarris

Blackman–Harris is a sum-of-cosines window, with the terms

\[W(x) = 0.35875 + 0.48829\cos(2 \pi x) + 0.14128\cos(4 \pi x) + 0.01168\cos(6 \pi x).\]

ElectricFields.BlackmanNuttallType
ElectricFields.BlackmanNuttall

Blackman–Nuttall is a sum-of-cosines window, with the terms

\[W(x) = 0.3635819 + 0.4891775\cos(2 \pi x) + 0.1365995\cos(4 \pi x) + 0.0106411\cos(6 \pi x).\]

ElectricFields.CascadedDispersiveElementType
CascadedDispersiveElement(elements)

Represents the combination of multiple DispersiveElements.

Example

julia> PhaseShift(6)*Chirp(austrip(5u"fs^2"), austrip(1.5u"eV"))
[PhaseShift(ϕ = 6.0000 rad) ∗ Chirp(b = 8545.5457 = 5.0000 fs², ω₀ = 0.0551 = 1.5000 eV)]
ElectricFields.ChirpType
Chirp(b, ω₀)

Represents chirp according to

\[H(\omega) = \exp[-\im b(ω-ω_0)^2].\]

Example

julia> Chirp(austrip(5u"fs^2"), austrip(1.5u"eV"))
Chirp(b = 8545.5457 = 5.0000 fs², ω₀ = 0.0551 = 1.5000 eV)
ElectricFields.ConstantFieldType
ConstantField(tmax, E₀)

The field amplitude of a constant field is defined as

\[F(t) = \begin{cases} E_0, & 0 \leq t \leq t_{\textrm{max}}, \\ 0, & \textrm{else}, \end{cases} \implies A(t) = \begin{cases} -E_0t, & 0 \leq t \leq t_{\textrm{max}}, \\ A(t_{\textrm{max}}), & t_{\textrm{max}} < t, \\ 0, & \textrm{else}. \end{cases}\]

Since the vector potential is non-zero at the end of the pulse, this is a non-propagating field, i.e. it does not correspond to a freely propagating pulse. It however corresponds to the field in an idealized capacitor, i.e. two plates of opposite charge.

Example

julia> @field(F) do
       I₀ = 1e13u"W/cm^2"
       tmax = 3.0u"fs"
       kind = :constant
       end
Constant field of
  - 124.0241 jiffies = 3.0000 fs duration, and
  - E₀ = 1.6880e-02 au = 8.6802 GV m^-1
ElectricFields.Cos²EnvelopeType
Cos²Envelope

\[f(t) = \begin{cases} \cos^2 \left(\frac{\pi t}{cT}\right), & -1 \leq t/(cT) \leq 1,\\ 0, & \textrm{else}, \end{cases}\]

where $c$ is the number of cycles from zero to zero of the $\cos^2$ envelope and $T$ the period time.

Required parameters:

  • λ|T|f|ν|ω|ħω,
  • cycles,
  • env=:cos² | env=:cos2.
ElectricFields.CrystalType
Crystal(material, d, R=I; ω₀)

Convenience constructor for Crystal; if a central angular frequency ω₀ is provided, the linear slope of the dispersion of the material is computed at ω₀. This slope is subsequently subtracted from the dispersion.

ElectricFields.CrystalType
Crystal(material, d, R[, ∂k∂ω₀ = 0])

Describes the dispersion through a general crystal of thickness d and orientation given by the rotation R, with different refractive indices, given by the components of material, along the different crystal axes.

In the special case where material has two components, we have a uniaxial crystal, which is birefringent, and material[1] is referred to as ordinary the refractive index, and material[2] as the extraordinary refractive index.

The optional ∂k∂ω₀ can be used to subtract the linear slope from the dispersion relation, i.e. remove the time shift, such that a pulse with central angular frequency ω₀ stays centred in the frame of reference.

ElectricFields.DelayedFieldType
DelayedField(a, t₀)

Represents a delayed copy of a, that appears at time t₀ instead of at 0, i.e. t₀>0 implies the field comes later.

Examples

julia> @field(A) do
           I₀ = 1.0
           T = 2.0
           σ = 3.0
           Tmax = 3.0
       end
Linearly polarized field with
  - I₀ = 1.0000e+00 au = 3.5094452e16 W cm^-2 =>
    - E₀ = 1.0000e+00 au = 514.2207 GV m^-1
    - A₀ = 0.3183 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0253 Ha = 689.2724 meV => α = 0.1013 Bohr = 5.3617 pm

julia> delay(A, 1u"fs")
Linearly polarized field with
  - I₀ = 1.0000e+00 au = 3.5094452e16 W cm^-2 =>
    - E₀ = 1.0000e+00 au = 514.2207 GV m^-1
    - A₀ = 0.3183 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0253 Ha = 689.2724 meV => α = 0.1013 Bohr = 5.3617 pm
  – delayed by 41.3414 jiffies = 1.0000 fs

julia> delay(A, 1.0)
Linearly polarized field with
  - I₀ = 1.0000e+00 au = 3.5094452e16 W cm^-2 =>
    - E₀ = 1.0000e+00 au = 514.2207 GV m^-1
    - A₀ = 0.3183 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0253 Ha = 689.2724 meV => α = 0.1013 Bohr = 5.3617 pm
  – delayed by 1.0000 jiffies = 24.1888 as

julia> delay(A, π*u"rad")
Linearly polarized field with
  - I₀ = 1.0000e+00 au = 3.5094452e16 W cm^-2 =>
    - E₀ = 1.0000e+00 au = 514.2207 GV m^-1
    - A₀ = 0.3183 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0253 Ha = 689.2724 meV => α = 0.1013 Bohr = 5.3617 pm
  – delayed by 1.0000 jiffies = 24.1888 as
ElectricFields.DiracCombType
DiracComb(frequencies)

Represents a Dirac frequency comb, where frequencies is a vector of Tuples: (ωᵢ,cᵢ), representing a frequency and an amplitude:

\[\mathrm{DC}(\omega) = \sum_i c_i \delta(\omega-\omega_i).\]

ElectricFields.DispersedFieldType
DispersedField(f, de, B, s)

Represents the field f dispersed through the DispersiveElement de. Since the dispersion is most efficiently calculated in the frequency domain, this is done via rfft/irfft, and the dispersed field is interpolated as the BSplineField B. The temporal span of the dispersed field is given by s.

ElectricFields.DispersedFieldMethod
DispersedField(f::DispersedField, de)

Stacking another DispersiveElement de onto an already dispersed field f is accomlished by multiplying the transfer functions together, with the earliest applied as the right-most factor:

\[H(\omega) = H_n(\omega) ... H_3(\omega)H_2(\omega)H_1(\omega).\]

ElectricFields.DispersedFieldMethod
DispersedField(f, de; spline_order=7, Bfs=20/period(F))

Convenience constructor for DispersedField that automatically finds the appropriate time span for f after it has been dispersed through de, and then fits a BSpline to the resultant vector potential, such that vector_potential and field_amplitude can be evaluated at arbitrary times. The number of knots can be influenced by the "sampling frequency" Bfs.

ElectricFields.EllipticalCarrierType
EllipticalCarrier

An elliptically polarized field (of which circularly and linearly polarized are special cases) is given in the canonical coordinate system as

\[\vec{C}(t) = \frac{1}{\sqrt{1+\xi^2}} \bmat{\xi\cos\theta\\0\\\sin\theta}, \quad \theta \defd \omega t + \phi,\]

i.e. with the principal axis of the ellipse along $z$.

ElectricFields.ExpKnotSetType
ExpKnotSet{k,ml,mr}(exponents, base, t, include0)

A knot set of order k and left and right multiplicities of ml and mr, respectively, whose knots are exponentially distributed according to t = base .^ exponents, optionally including 0 as the left endpoint.

ElectricFields.ExpKnotSetMethod
ExpKnotSet(k, a, b, N[, ml=k, mr=k, base=10, include0=true])

Construct an order-k knot spanning from base^a to base^b in N intervals, optionally including 0 as the left endpoint.

Examples

julia> ElectricFields.ExpKnotSet(2, -4, 2, 7)
10-element ElectricFields.ExpKnotSet{2, 2, 2, Float64, StepRangeLen{Float64, Base.TwicePrecision{Float64}, Base.TwicePrecision{Float64}, Int64}, Vector{Float64}}:
   0.0
   0.0
   0.0001
   0.001
   0.010000000000000002
   0.1
   1.0
  10.0
 100.0
 100.0
ElectricFields.FixedCarrierType
FixedCarrier <: LinearCarrier

The carrier is fixed in the sense that the instantaneous frequency is constant throughout the pulse, i.e. no chirp/dispersion:

\[\Im\{\exp[\im(\omega t + \phi)]\} = \sin(\omega t + \phi).\]

ElectricFields.GaussianEnvelopeType
GaussianEnvelope

A Gaussian pulse is given by

\[f(t) = \exp\left(-\frac{t^2}{2\sigma^2}\right),\]

where the standard deviation $σ$ is related to the FWHM duration τ of the intensity envelope as

\[\sigma = \frac{\tau}{2\sqrt{2\ln 2}}.\]

Since we define all fields in terms of the vector potential $A(t)$, and the instantaneous_intensity is given by $\abs{-\partial_t A(t)}^2$, we have to find a coefficient $α$ such that the field amplitude

\[F(t) \sim -\partial_t \exp(-\alpha t^2) \sin(\omega t + \phi)\]

has an intensity envelope with the desired FWHM; we do this iteratively in gaussian_common!. This is mainly important for ultrashort pulses, since for long pulse durations we can approximate $\exp(-\alpha t^2) \sim 1$ such that only the carrier contributes to the time derivative.

Since a Gaussian never ends, we specify how many $σ$ we require; the resulting time window will be rounded up to an integer amount of cycles of the fundamental.

Required parameters:

  • λ|T|f|ν|ω|ħω,
  • τ|σ|α,
  • σmax|tmax|Tmax,
  • env=:gauss (optional, since GaussianEnvelope is the default envelope).
ElectricFields.HammingType
ElectricFields.Hamming

Hamming is a sum-of-cosines window, with the terms

\[W(x) = 0.5434782608695652 + 0.45652173913043476\cos(2 \pi x).\]

ElectricFields.HannType
ElectricFields.Hann

Hann is a sum-of-cosines window, with the terms

\[W(x) = 0.5 + 0.5\cos(2 \pi x).\]

ElectricFields.IsotropicMediumType
IsotropicMedium(material, d[, ∂k∂ω₀ = 0])

Describes the dispersion through an isotropic medium of thickness d (such as an atomic gas), where the refractive index, given by material, is the same in all directions.

The optional ∂k∂ω₀ can be used to subtract the linear slope from the dispersion relation, i.e. remove the time shift, such that a pulse with central angular frequency ω₀ stays centred in the frame of reference.

ElectricFields.IsotropicMediumMethod
IsotropicMedium(material, d; ω₀)

Convenience constructor for IsotropicMedium; if a central angular frequency ω₀ is provided, the linear slope of the dispersion of the material is computed at ω₀. This slope is subsequently subtracted from the dispersion.

ElectricFields.KaiserType
Kaiser(α)

\[W(x) = \frac{I_0\left[\pi\alpha\sqrt{1 - (2x)^2}\right]}{I_0(\pi\alpha)},\]

where $I_0$ is the 0th order modified Bessel function of the first kind.

ElectricFields.LinearCarrierType
LinearCarrier <: AbstractCarrier

Any carrier which is one-dimensional, i.e. linearly polarized (by convention along $z$).

ElectricFields.LinearFieldType
LinearField

Linearly polarized field, i.e. the field amplitude is scalar. Consisting of a peak vector potential, a vector envelope, and a carrier, which has to be a LinearCarrier.

Examples

julia> @field(A) do
       I₀ = 1e14u"W/cm^2"
       λ = 800.0u"nm"
       τ = 6.2u"fs"
       tmax = 20.0u"fs"
       end
Linearly polarized field with
  - I₀ = 2.8495e-03 au = 1.0e14 W cm^-2 =>
    - E₀ = 5.3380e-02 au = 27.4492 GV m^-1
    - A₀ = 0.9372 au
  – a Fixed carrier @ λ = 800.0000 nm (T = 2.6685 fs, ω = 0.0570 Ha = 1.5498 eV, f = 374.7406 THz)
  – and a Gaussian envelope of duration 6.2000 fs (intensity FWHM; ±8.11σ)
  – and a bandwidth of 0.0108 Ha = 294.3469 meV ⟺ 71.1728 THz ⟺ 2871.2568 Bohr = 151.9404 nm
  – Uₚ = 0.2196 Ha = 5.9759 eV => α = 16.4562 Bohr = 870.8242 pm

julia> @field(B) do
       I₀ = 0.05
       ω = 1.0
       ramp = 1.0
       flat = 3.0
       env = :trapezoidal
       end
Linearly polarized field with
  - I₀ = 5.0000e-02 au = 1.7547226e15 W cm^-2 =>
    - E₀ = 2.2361e-01 au = 114.9832 GV m^-1
    - A₀ = 0.2236 au
  – a Fixed carrier @ λ = 45.5634 nm (T = 151.9830 as, ω = 1.0000 Ha = 27.2114 eV, f = 6.5797 PHz)
  – and a /1‾3‾1\ cycles trapezoidal envelope, with linear ramps
  – and a bandwidth of Inf Ha = Inf eV ⟺ Inf Hz ⟺ Inf Bohr = Inf m
  – Uₚ = 0.0125 Ha = 340.1423 meV => α = 0.2236 Bohr = 11.8328 pm
ElectricFields.LinearKnotSetType
LinearKnotSet(k, a, b, N[, ml=k, mr=k])

Construct an order-k linear knot set spanning from a to b, with N intervals.

Examples

julia> ElectricFields.LinearKnotSet(2, 0, 1, 3)
6-element ElectricFields.LinearKnotSet{2, 2, 2, Float64, StepRangeLen{Float64, Base.TwicePrecision{Float64}, Base.TwicePrecision{Float64}, Int64}}:
 0.0
 0.0
 0.3333333333333333
 0.6666666666666666
 1.0
 1.0
ElectricFields.LinearKnotSetType
LinearKnotSet{k,ml,mr}(t)

A knot set of order k and left and right multiplicities of ml and mr, respectively, whose knots are uniformly distributed according to the range t. The interior basis functions are thus Cᵏ⁻²-continuous.

ElectricFields.LinearTransverseFieldType

LinearTransverseField

Linearly polarized field, but explicitly in 3d, i.e the polarization vector is confined to a line in the plane that is perpendicular to the direction of propagation. If the propagation vector is not rotated, this plane is the $z-x$, with $z$ by convention being the principal polarization axis, and $y$ is the direction of propagation.

This is mostly used as a proxy, if one wants to calculate using a manifestly linearly polarized field in 3d. See also LinearTransverseCarrier.

ElectricFields.NegatedFieldType
NegatedField(a)

Represents a field whose vector_potential is the negative of that of a.

Example

julia> @field(A) do
           I₀ = 1.0
           T = 2.0
           σ = 3.0
           Tmax = 3.0
       end
Linearly polarized field with
  - I₀ = 1.0000e+00 au = 3.5094452e16 W cm^-2 =>
    - E₀ = 1.0000e+00 au = 514.2207 GV m^-1
    - A₀ = 0.3183 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0253 Ha = 689.2724 meV => α = 0.1013 Bohr = 5.3617 pm

julia> B = -A
ElectricFields.NegatedField{ElectricFields.LinearField{ElectricFields.FixedCarrier{Quantity{Float64, 𝐋, Unitful.FreeUnits{(Eₕ^-1, ħ, c), 𝐋, nothing}}, Quantity{Float64, 𝐓, Unitful.FreeUnits{(Eₕ^-1, ħ), 𝐓, nothing}}, Float64, Int64}, ElectricFields.GaussianEnvelope{Float64}, Float64}}(Linearly polarized field with
  - I₀ = 1.0000e+00 au = 3.5094452e16 W cm^-2 =>
    - E₀ = 1.0000e+00 au = 514.2207 GV m^-1
    - A₀ = 0.3183 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0253 Ha = 689.2724 meV => α = 0.1013 Bohr = 5.3617 pm)

julia> field_amplitude(A, 0.5)
0.008830294133325867

julia> field_amplitude(B, 0.5)
-0.008830294133325867

julia> field_amplitude(A-A, 0.5)
0.0
ElectricFields.NuttallType
ElectricFields.Nuttall

Nuttall is a sum-of-cosines window, with the terms

\[W(x) = 0.355768 + 0.487396\cos(2 \pi x) + 0.144232\cos(4 \pi x) + 0.012604\cos(6 \pi x).\]

ElectricFields.PaddedFieldType
  PaddedField(field, a, b)

Wrapper around any electric field, padded with a units of time before, and b units of time after the ordinary span of the field.

Example

julia> @field(A) do
           I₀ = 1.0
           T = 2.0
           σ = 3.0
           Tmax = 3.0
       end
Linearly polarized field with
  - I₀ = 1.0000e+00 au = 3.5094452e16 W cm^-2 =>
    - E₀ = 1.0000e+00 au = 514.2207 GV m^-1
    - A₀ = 0.3183 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0253 Ha = 689.2724 meV => α = 0.1013 Bohr = 5.3617 pm

julia> B = PaddedField(A, 10.0, 30.0)
Padding before 10.0000 jiffies = 241.8884 as and after 30.0000 jiffies = 725.6653 as of
Linearly polarized field with
  - I₀ = 1.0000e+00 au = 3.5094452e16 W cm^-2 =>
    - E₀ = 1.0000e+00 au = 514.2207 GV m^-1
    - A₀ = 0.3183 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0253 Ha = 689.2724 meV => α = 0.1013 Bohr = 5.3617 pm

julia> span(A), span(B)
(-6.0 .. 6.0, -16.0 .. 36.0)
ElectricFields.PhaseShiftType
PhaseShift(ϕ)

Represents a phase shift according to

\[H(\omega) = \exp(-\im\phi).\]

Example

julia> PhaseShift(6)
PhaseShift(ϕ = 6.0000 rad)
ElectricFields.RampType
Ramp(tmax, E₀, f, name)

The field amplitude of a ramp is defined as

\[F(t) = \begin{cases} E_0f'(\tau), & 0 \leq t_{\textrm{max}}, \\ 0, \textrm{else}, \end{cases} \implies A(t) = \begin{cases} -E_0t_{\textrm{max}}f(\tau), & 0 \leq t \leq t_{\textrm{max}},\\ A(t_{\textrm{max}}), & t_{\textrm{max}} < t,\\ 0, & \textrm{else}, \end{cases} \quad \tau = \frac{t}{t_{\textrm{max}}}.\]

To define a new ramp, one thus only needs to define one function on the unit interval, $f(\tau)$ whose derivative $f'(\tau)$ rises from $0$ to $1$. Similar to ConstantField, Ramp is a non-propagating field, but is realizable in e.g. a capacitor.

Three kinds of ramps are predefined, :linear_ramp, :parabolic_ramp, :sin²_ramp (with the alias :sin2_ramp).

Examples

julia> @field(F) do
           E₀ = 1.0
           tmax = 4.0u"fs"
           kind = :sin²_ramp
       end
sin² up-ramp of
  - 165.3655 jiffies = 4.0000 fs duration, and
  - E₀ = 1.0000e+00 au = 514.2207 GV m^-1

julia> @field(F) do
           E₀ = 1.0
           tmax = 4.0u"fs"
           kind = :linear_ramp
           ramp = :down
       end
Linear down-ramp of
  - 165.3655 jiffies = 4.0000 fs duration, and
  - E₀ = 1.0000e+00 au = 514.2207 GV m^-1
ElectricFields.RectType
Rect

Rectangular window:

\[W(x) = \begin{cases} 1, & |x| \le \frac{1}{2}, \\ 0, & \textrm{else}. \end{cases}\]

ElectricFields.SellmeierType
Sellmeier(A, B, p, C, D, q)

The Sellmeier equations are used to describe dispersion in glass using a series of resonances, here generalized to support variations found in the literature:

\[n^2(\lambda) = 1 + A + \sum_i \frac {B_i\lambda^{p_i}}{\lambda^2-C_i} + \sum_j D_j \lambda^{q_j}.\]

ElectricFields.SellmeierMethod
Sellmeier(B, C)

Convenience constructor for Sellmeier for the canonical form of the Sellmeier equation:

\[n^2(\lambda) = 1 + \sum_i \frac {B_i\lambda^2}{\lambda^2-C_i}.\]

ElectricFields.SumFieldType
SumField(a, b)

The linear combination of two fields a and b.

Example

julia> @field(A) do
           I₀ = 1.0
           T = 2.0
           σ = 3.0
           Tmax = 3.0
       end
Linearly polarized field with
  - I₀ = 1.0000e+00 au = 3.5094452e16 W cm^-2 =>
    - E₀ = 1.0000e+00 au = 514.2207 GV m^-1
    - A₀ = 0.3183 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0253 Ha = 689.2724 meV => α = 0.1013 Bohr = 5.3617 pm

julia> @field(B) do
           I₀ = 0.5
           T = 1.0
           σ = 3.0
           Tmax = 3.0
       end
Linearly polarized field with
  - I₀ = 5.0000e-01 au = 1.7547226e16 W cm^-2 =>
    - E₀ = 7.0711e-01 au = 363.6089 GV m^-1
    - A₀ = 0.1125 au
  – a Fixed carrier @ λ = 7.2516 nm (T = 24.1888 as, ω = 6.2832 Ha = 170.9742 eV, f = 41.3414 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±1.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 8.5598 Bohr = 452.9627 pm
  – Uₚ = 0.0032 Ha = 86.1591 meV => α = 0.0179 Bohr = 947.8211 fm

julia> A+B
┌ Linearly polarized field with
│   - I₀ = 1.0000e+00 au = 3.5094452e16 W cm^-2 =>
│     - E₀ = 1.0000e+00 au = 514.2207 GV m^-1
│     - A₀ = 0.3183 au
│   – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
│   – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
│   – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
│   – Uₚ = 0.0253 Ha = 689.2724 meV => α = 0.1013 Bohr = 5.3617 pm
⊕
│ Linearly polarized field with
│   - I₀ = 5.0000e-01 au = 1.7547226e16 W cm^-2 =>
│     - E₀ = 7.0711e-01 au = 363.6089 GV m^-1
│     - A₀ = 0.1125 au
│   – a Fixed carrier @ λ = 7.2516 nm (T = 24.1888 as, ω = 6.2832 Ha = 170.9742 eV, f = 41.3414 PHz)
│   – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±1.00σ)
│   – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 8.5598 Bohr = 452.9627 pm
└   – Uₚ = 0.0032 Ha = 86.1591 meV => α = 0.0179 Bohr = 947.8211 fm
ElectricFields.TransverseCarrierType
TransverseCarrier <: AbstractCarrier

Any carrier which is two-dimensional, i.e. polarized in the plane perpendicular to the direction of propagation.

ElectricFields.TransverseFieldType
TransverseField

Transversely polarized field, i.e the polarization vector is confined to the plane that is perpendicular to the direction of propagation. If the propagation vector is not rotated, this plane is the $z-x$, with $z$ by convention being the principal polarization axis, and $y$ is the direction of propagation.

The formal definition is

\[\vec{A}(t)\defd A_0 f(t) \Im\{\vec{J}\exp[\im(\omega t + \phi)]\},\]

using the complex-valued Jones vector $\vec{J}$, but in practice, we use a rotation matrix $\mat{R}$ and the Carriers give the amplitudes in the canonical, unrotated coordinate sytem:

\[\vec{A}(t) = A_0 f(t) \mat{R}\vec{C}(t).\]

Examples

julia> # Circularly polarized field

julia> @field(A) do
       I₀ = 1e14u"W/cm^2"
       λ = 800.0u"nm"
       τ = 6.2u"fs"
       tmax = 20.0u"fs"
       ξ = 1.0
       end
Transversely polarized field with
  - I₀ = 2.8495e-03 au = 1.0e14 W cm^-2 =>
    - E₀ = 5.3380e-02 au = 27.4492 GV m^-1
    - A₀ = 0.9372 au
  – a Elliptical carrier with ξ = 1.00 (RCP) @ λ = 800.0000 nm (T = 2.6685 fs, ω = 0.0570 Ha = 1.5498 eV, f = 374.7406 THz)
  – and a Gaussian envelope of duration 6.2000 fs (intensity FWHM; ±8.11σ)
  – and a bandwidth of 0.0108 Ha = 294.3469 meV ⟺ 71.1728 THz ⟺ 2871.2568 Bohr = 151.9404 nm
  – Uₚ = 0.2196 Ha = 5.9759 eV => α = 16.4562 Bohr = 870.8242 pm

julia> # Linearly polarized field, but explicitly in 3D

julia> @field(B) do
       I₀ = 1e14u"W/cm^2"
       λ = 800.0u"nm"
       τ = 6.2u"fs"
       tmax = 20.0u"fs"
       kind = :transverse
       end
Transversely polarized field with
  - I₀ = 2.8495e-03 au = 1.0e14 W cm^-2 =>
    - E₀ = 5.3380e-02 au = 27.4492 GV m^-1
    - A₀ = 0.9372 au
  – a LinearTransverseCarrier: Fixed carrier @ λ = 800.0000 nm (T = 2.6685 fs, ω = 0.0570 Ha = 1.5498 eV, f = 374.7406 THz)
  – and a Gaussian envelope of duration 6.2000 fs (intensity FWHM; ±8.11σ)
  – and a bandwidth of 0.0108 Ha = 294.3469 meV ⟺ 71.1728 THz ⟺ 2871.2568 Bohr = 151.9404 nm
  – Uₚ = 0.2196 Ha = 5.9759 eV => α = 16.4562 Bohr = 870.8242 pm

julia> # Linearly polarized field, rotated

julia> @field(C) do
       I₀ = 1e14u"W/cm^2"
       λ = 800.0u"nm"
       τ = 6.2u"fs"
       tmax = 20.0u"fs"
       rotation = π/3, [0,0,1]
       end
Transversely polarized field with
  - I₀ = 2.8495e-03 au = 1.0e14 W cm^-2 =>
    - E₀ = 5.3380e-02 au = 27.4492 GV m^-1
    - A₀ = 0.9372 au
  – a LinearTransverseCarrier: Fixed carrier @ λ = 800.0000 nm (T = 2.6685 fs, ω = 0.0570 Ha = 1.5498 eV, f = 374.7406 THz)
  – a Gaussian envelope of duration 6.2000 fs (intensity FWHM; ±8.11σ)
  – and a rotation of 0.33π about [0.000, 0.000, 1.000]
  – and a bandwidth of 0.0108 Ha = 294.3469 meV ⟺ 71.1728 THz ⟺ 2871.2568 Bohr = 151.9404 nm
  – Uₚ = 0.2196 Ha = 5.9759 eV => α = 16.4562 Bohr = 870.8242 pm
ElectricFields.TrapezoidalEnvelopeType
TrapezoidalEnvelope

This is a very simple piecewise linear function:

\[f(t)= \begin{cases} r/r_{\textrm{up}}, & 0 \leq r < r_{\textrm{up}},\\ 1, & r_{\textrm{up}} \leq r < r_{\textrm{up}} + r_{\textrm{flat}}, \\ 1 - \frac{r-r_{\textrm{up}}-r_{\textrm{flat}}}{r_{\textrm{down}}}, & r_{\textrm{up}} + r_{\textrm{flat}} \leq r < r_{\textrm{up}} + r_{\textrm{flat}} + r_{\textrm{down}}, \\ 0, & \textrm{else}, \end{cases} \quad r \defd t/T,\]

where $T$ is the period time of the carrier.

Required parameters:

  • λ|T|f|ν|ω|ħω,
  • ramp | (ramp_up & ramp_down),
  • flat,
  • env=:trapezoidal | env=:tophat.

Beware that this envelope can introduce artifacts at the ends of the pulse, such that the electric field is non-vanishing, depending of e.g. the phase of the carrier.

The shape of the ramp can be influenced by chaining $g[f(t)]$, where $g(x)$ should map smoothly to $[0,1]$ for $x\in[0,1]$. Currently implemented ramps are

  • ramp_kind = :linear $\implies g(x) = x$ (default)
  • ramp_kind = :sin² | ramp_kind = :sin2 $\implies g(x) = \sin^2(\pi x/2)$
ElectricFields.TruncatedGaussianEnvelopeType
TruncatedGaussianEnvelope

Since a Gaussian function never ends, this envelope adds a soft truncation over a time interval:

\[f(t)= \begin{cases} \exp(-\alpha t^2), & \abs{t} \leq t_{\textrm{off}},\\ \exp\left\{ -\alpha\left[ t_{\textrm{off}} + \frac{2}{\pi}(t_{\textrm{max}} - t_{\textrm{off}}) \tan\left( \frac{\pi}{2} \frac{\abs{t} - t_{\textrm{off}}}{t_{\textrm{max}} - t_{\textrm{off}}} \right) \right]^2 \right\}, & t_{\textrm{off}} < \abs{t} \leq t_{\textrm{max}}, \\ 0, & t_{\textrm{max}} < \abs{t}. \end{cases}\]

This is Eq. (72) of

  • Patchkovskii, S., & Muller, H. (2016). Simple, accurate, and efficient implementation of 1-electron atomic time-dependent Schrödinger equation in spherical coordinates. Computer Physics Communications, 199, 153–169. 10.1016/j.cpc.2015.10.014

Required parameters:

  • λ|T|f|ν|ω|ħω,
  • τ|σ|α,
  • toff,
  • σmax|tmax|Tmax,

env=:trunc_gauss.

Beyond this, everything else is the same as for GaussianEnvelope.

ElectricFields.UnitVectorType
UnitVector{T}(N, k)

Helper vector type of length N where the kth element is one(T) and all the others zero(T).

ElectricFields.WindowedFieldType
WindowedField(field, a, b)

Wrapper around any electric field, windowed such that it is zero outside the time interval a..b.

Example

julia> @field(A) do
           I₀ = 1.0
           T = 2.0
           σ = 3.0
           Tmax = 3.0
       end
Linearly polarized field with
  - I₀ = 1.0000e+00 au = 3.5094452e16 W cm^-2 =>
    - E₀ = 1.0000e+00 au = 514.2207 GV m^-1
    - A₀ = 0.3183 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0253 Ha = 689.2724 meV => α = 0.1013 Bohr = 5.3617 pm

julia> B = WindowedField(A, -3, 5)
Window from -3.0000 jiffies = -72.5665 as to 5.0000 jiffies = 120.9442 as of
Linearly polarized field with
  - I₀ = 1.0000e+00 au = 3.5094452e16 W cm^-2 =>
    - E₀ = 1.0000e+00 au = 514.2207 GV m^-1
    - A₀ = 0.3183 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0253 Ha = 689.2724 meV => α = 0.1013 Bohr = 5.3617 pm

julia> span(A), span(B)
(-6.0 .. 6.0, -3.0 .. 5.0)

julia> field_amplitude(A, -4)
-0.6395632362683398

julia> field_amplitude(B, -4)
0.0
AbstractFFTs.fftMethod
fft(f::AbstractField, t::AbstractRange)

Compute the FFT of the field_amplitude of the field f sampled on the uniform temporal grid t.

AbstractFFTs.irfftMethod
irfft(F, t)

Compute the IRFFT of the spectral amplitudes F, with the end-result resolved on the time axis t, in the case of real-valued signals.

AbstractFFTs.rfftMethod
rfft(f::AbstractField, t::AbstractRange)

Compute the RFFT of the field_amplitude of the field f sampled on the uniform temporal grid t, in the case of real-valued signals.

Base.firstMethod
first(t)

Return the first knot of t.

Base.getindexMethod
getindex(t, i)

Return the ith knot of the knot set t, accounting for the multiplicities of the endpoints.

Examples

julia> ElectricFields.LinearKnotSet(3, 0, 1, 3)[2]
0.0

julia> ElectricFields.LinearKnotSet(3, 0, 1, 3, 1, 1)[2]
0.3333333333333333
Base.lastMethod
last(t)

Return the last knot of t.

Base.lengthMethod
length(t)

Return the number of knots of t.

Examples

julia> length(ElectricFields.LinearKnotSet(3, 0, 1, 3))
8

julia> length(ElectricFields.LinearKnotSet(3, 0, 1, 3, 1, 1))
4
ElectricFields.assert_multiplicitiesMethod
assert_multiplicities(k,ml,mr,t)

Assert that the multiplicities at the endpoints, ml and mr, respectively, are consistent with the order k. Also check that the amount of knots in the knot set t are enough to support the requested order k.

ElectricFields.calc_params!Method
calc_params!(field_params)

This function performs the calculation of different quantities from the information provided.

ElectricFields.chirpFunction
chirp(f, b, ω₀=photon_energy(f))

Returns the field resulting from applying a Chirp to the field f.

Example

julia> @field(F) do
           I₀ = 1.0
           T = 2.0u"fs"
           σ = 3.0
           Tmax = 3.0
       end
Linearly polarized field with
  - I₀ = 1.0000e+00 au = 3.5094452e16 W cm⁻² =>
    - E₀ = 1.0000e+00 au = 514.2207 GV m⁻¹
    - A₀ = 13.1594 au
  – a Fixed carrier @ λ = 599.5849 nm (T = 2.0000 fs, ω = 0.0760 Ha = 2.0678 eV, f = 500.0000 THz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±82.68σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 58518.2144 Bohr = 3.0967 μm
  – Uₚ = 43.2922 Ha = 1.1780 keV => α = 173.1690 Bohr = 9.1637 nm

julia> chirp(F, austrip(1u"fs^2"))
DispersedField:
Linearly polarized field with
  - I₀ = 1.0000e+00 au = 3.5094452e16 W cm⁻² =>
    - E₀ = 1.0000e+00 au = 514.2207 GV m⁻¹
    - A₀ = 13.1594 au
  – a Fixed carrier @ λ = 599.5849 nm (T = 2.0000 fs, ω = 0.0760 Ha = 2.0678 eV, f = 500.0000 THz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±82.68σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 58518.2144 Bohr = 3.0967 μm
  – Uₚ = 43.2922 Ha = 1.1780 keV => α = 173.1690 Bohr = 9.1637 nm
  – dispersed through Chirp(b = 1709.1091 = 1.0000 fs², ω₀ = 0.0760 = 2.0678 eV)
ElectricFields.convolutionMethod
convolution(f̂::Function, dc::DiracComb, ω)

Evaluate the convolution between the function (assumed to be the Fourier transform of a function f) with the DiracComb dc. This is used to implement the Fourier transform of a function product $f(t)g(t)$, where $g(t)$ is a sum of monochromatic waves:

\[f(t)g(t) \rightsquigarrow \frac{1}{\sqrt{2\pi}} (\hat{f}\star\hat{g})(\omega)\]

ElectricFields.deBoorFunction
deBoor(t, c, x[, i[, m=0]])

Evaluate the spline given by the knot set t and the set of control points c at x using de Boor's algorithm. i is the index of the knot interval containing x. If m≠0, calculate the mth derivative at x instead.

ElectricFields.dimensionsFunction
dimensions(f)

Return the number of dimensions of the field f. See also polarization.

Examples

julia> @field(F) do
           I₀ = 2.0
           T = 2.0
           σ = 3.0
           Tmax = 3.0
       end
Linearly polarized field with
  - I₀ = 2.0000e+00 au = 7.0188904e16 W cm^-2 =>
    - E₀ = 1.4142e+00 au = 727.2178 GV m^-1
    - A₀ = 0.4502 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0507 Ha = 1.3785 eV => α = 0.1433 Bohr = 7.5826 pm

julia> dimensions(F)
1

julia> @field(F) do
           I₀ = 2.0
           T = 2.0
           σ = 3.0
           Tmax = 3.0
           ξ = 1.0
       end
Transversely polarized field with
  - I₀ = 2.0000e+00 au = 7.0188904e16 W cm^-2 =>
    - E₀ = 1.4142e+00 au = 727.2178 GV m^-1
    - A₀ = 0.4502 au
  – a Elliptical carrier with ξ = 1.00 (RCP) @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0507 Ha = 1.3785 eV => α = 0.1433 Bohr = 7.5826 pm

julia> dimensions(F)
3
ElectricFields.fftωMethod
fftω(t::AbstractRange)

Return the angular frequency grid corresponding to uniform sampling in time with the tempolar grid t.

ElectricFields.field_amplitude_spectrumMethod
field_amplitude_spectrum(f::AbstractField, ω)

Compute the analytic Fourier transform of the field_amplitude of the field f at the angular frequency ω, using the Fourier identity

\[\vec{F}(t) = -\partial_t \vec{A}(t) \iff \hat{\vec{F}}(\omega) = -\im\omega \hat{\vec{A}}(\omega)\]

ElectricFields.find_intervalMethod
find_interval(t, x[, i=ml])

Find the interval in the knot set t that includes x, starting from interval i (which by default is the first non-zero interval of the knot set). The search complexity is linear, but by storing the result and using it as starting point for the next call to find_interval, the knot set need only be traversed once.

ElectricFields.find_time_spanMethod
find_time_span(f, de[, fs]; max_iter=10, ξ=2.0, tol)

Find a suitable time span $[a,b]$ such that when f is dispersed through the DispersiveElement de, the compact support of the resulting field is contained in the time interval. This is done by successively multiplying the time span on which F is evaluated ξ before the RFFT, until the IRFFT has converged.

ElectricFields.fluenceMethod
fluence(f)

Compute the fluence of the field f, i.e.

\[\frac{1}{\hbar\omega} \int\diff{t} I(t),\]

where $I(t)$ is the intensity envelope of the pulse.

ElectricFields.free_oscillation_amplitudeMethod
free_oscillation_amplitude(F)

Compute the free oscillation amplitude of an electric field F, i.e. the mean excursion length during one cycle of the field, defined as

\[\alpha \defd \frac{F}{\omega^2}\]

where F is the peak amplitude, i.e. this is defined for one cycle of a monochrome field.

ElectricFields.gaussian_common!Method
gaussian_common!(field_params, carrier[; Tmax_rounder, verbosity])

Compute parameters common to Gaussian envelopes, i.e. GaussianEnvelope and TruncatedGaussianEnvelope; most importantly, given an intensity FWHM or σ, we need to figure out the coefficient $\alpha$ for the envelope of the vector potential $\exp(-\alpha t^2)$, such that the electric field amplitude and intensity have the desired durations. The optional function Tmax_rounder determines if the envelop should be extended to encompass an integer amount of cycles of the carrier (default).

ElectricFields.keldyshMethod
keldysh(f, Iₚ)

The Keldysh parameter relates the strength of a dynamic electric field to that of the binding potential of an atom. It is given by

\[\gamma = \sqrt{\frac{I_p}{2U_p}},\]

where $I_p$ is the ionization potential of the atom and $U_p$ is the ponderomotive potential of the dynamic field.

ElectricFields.nfftMethod
nfft(y, ts...)

Normalized FFT of y with respect to the axes ts. We use the symmetric normalization $(2\pi)^{-N/2}$ where $N$ is the number of dimensions.

ElectricFields.nfftMethod
nfft(f::AbstractField, t::AbstractRange)

Compute the symmetrically normalized FFT of the field_amplitude of the field f sampled on the uniform temporal grid t.

ElectricFields.nonempty_intervalsMethod
nonempty_intervals(t)

Return the indices of all intervals of the knot set t that are non-empty.

Examples

julia> ElectricFields.nonempty_intervals(ElectricFields.ArbitraryKnotSet(3, [0.0, 1, 1, 3, 4, 6], 1, 3))
4-element Vector{Int64}:
 1
 3
 4
 5
ElectricFields.numfunctionsMethod
numfunctions(t)

Returns the number of basis functions generated by knot set t.

Examples

julia> ElectricFields.numfunctions(ElectricFields.LinearKnotSet(3, 0, 1, 2))
4

julia> ElectricFields.numfunctions(ElectricFields.LinearKnotSet(5, 0, 1, 2))
6
ElectricFields.numintervalsMethod
numinterval(t)

Returns the number of intervals generated by the knot set t.

Examples

julia> ElectricFields.numintervals(ElectricFields.LinearKnotSet(3, 0, 1, 2))
2
ElectricFields.polarizationFunction
polarization(f)

Return the polarization kind of the field f. See also dimensions.

Examples

julia> @field(F) do
           I₀ = 2.0
           T = 2.0
           σ = 3.0
           Tmax = 3.0
       end
Linearly polarized field with
  - I₀ = 2.0000e+00 au = 7.0188904e16 W cm^-2 =>
    - E₀ = 1.4142e+00 au = 727.2178 GV m^-1
    - A₀ = 0.4502 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0507 Ha = 1.3785 eV => α = 0.1433 Bohr = 7.5826 pm

julia> polarization(F)
LinearPolarization()

julia> @field(F) do
           I₀ = 2.0
           T = 2.0
           σ = 3.0
           Tmax = 3.0
           ξ = 1.0
       end
Transversely polarized field with
  - I₀ = 2.0000e+00 au = 7.0188904e16 W cm^-2 =>
    - E₀ = 1.4142e+00 au = 727.2178 GV m^-1
    - A₀ = 0.4502 au
  – a Elliptical carrier with ξ = 1.00 (RCP) @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV, f = 20.6707 PHz)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
  – and a bandwidth of 0.3925 Ha = 10.6797 eV ⟺ 2.5823 PHz ⟺ 34.2390 Bohr = 1.8119 nm
  – Uₚ = 0.0507 Ha = 1.3785 eV => α = 0.1433 Bohr = 7.5826 pm

julia> polarization(F)
ArbitraryPolarization()
ElectricFields.ponderomotive_potentialMethod
ponderomotive_potential(f)

Return the ponderomotive potential $U_p$, which is the cycle-average quiver energy of a free electron in an electromagnetic field f. It is given by

\[U_p = \frac{e^2E_0^2}{4m\omega^2}=\frac{2e^2}{c\varepsilon_0m}\times\frac{I}{4\omega^2},\]

or, in atomic units,

\[U_p = \frac{I}{4\omega^2}.\]

ElectricFields.rfftωMethod
rfftω(t::AbstractRange)

Return the angular frequency grid corresponding to uniform sampling in time with the tempolar grid t, in the case of real-valued signals.

ElectricFields.spanMethod
span(f::ApodizedField)

The span of f is given by the intersection of span(parent(f)) and the apodizing interval a .. b.

ElectricFields.spectrumMethod
spectrum(env::AbstractGaussianEnvelope)

Gaussians belong to the Schwartz class, i.e. functions who, under Fourier transform, are mapped back to the same space. That is to say, the Fourier transform of a Gaussian is a Gaussian:

\[\exp(-\alpha t^2) \leftrightarrow \frac{1}{\sqrt{2\alpha}} \exp\left(-\frac{\omega^2}{4\alpha}\right).\]

Comparing with the above, we find that the spectral standard deviation

\[\Omega = \sqrt{2\alpha} = \frac{2\sqrt{\ln 2}}{\tau},\]

and the Gaussian function in the spectral domain is thus

\[E(\omega) = \frac{E_0\tau}{2\sqrt{\ln 2}} \exp\left[-\frac{(\omega\tau)^2}{8\ln2}\right].\]

ElectricFields.stepsFunction
steps(f[, fs])

Return number of time steps for the field f, optionally supplying a sampling frequency.

We construct the time axis such that the highest frequency component of the field is resolved. By the Nyquist sampling theorem, we need minimum $f_s>2f_{\textrm{max}}$, but to be on the safe side, we use, as default, $f_s=100f_{\textrm{max}}$. This also makes plots nicer.

If fs<:Integer, the sampling frequency is computed as fs/T, i.e. fs steps per cycle.

ElectricFields.timeaxisFunction
timeaxis(f[, fs])

Construct a time axis for the field f, covering the span of the envelope in the number of steps given by the sample frequency fs.

ElectricFields.timeaxisMethod
timeaxis(f, (N,dt))

Alternate interface that constructs the time axis with N steps starting from first(span(f)), with intervals of dt.

ElectricFields.within_supportMethod
within_support(x, t, j)

Return the indices of the elements of x that lie withing the compact support of the jth basis function (enumerated 1..n), given the knot set t. For each index of x that is covered, the index k of the interval within which x[i] falls is also returned.

ElectricFields.@fieldMacro
@field(name) do ... end

Frontend macro for creating an electric field and storing it in the variable name.

Example

julia> @field(F) do
           I₀ = 2.0
           T = 2.0
           σ = 3.0
           Tmax = 3.0
       end
Linearly polarized field with
  - I₀ = 2.0000e+00 au = 7.0188904e16 W cm⁻² =>
    - E₀ = 1.4142e+00 au = 727.2178 GV m⁻¹
    - A₀ = 0.4502 au
  – a Fixed carrier @ λ = 14.5033 nm (T = 48.3777 as, ω = 3.1416 Ha = 85.4871 eV)
  – and a Gaussian envelope of duration 170.8811 as (intensity FWHM; ±2.00σ)
ElectricFields.@namespace!Macro
@namespace!(exprs, params)

This macro uses the dictionary params as a "namespace", i.e. all symbols are assumed to be keys in this dictionary. We need to escape the generated expression tree to modify in the scope of the caller, not the global scope.