Public Interface
Documentation of LongwaveModePropagator.jl's exported structs and functions.
Contents
Basic Functions
LongwaveModePropagator.propagate — Functionpropagate(waveguide::HomogeneousWaveguide, tx::Emitter, rx::AbstractSampler;
modes::Union{Nothing,Vector{EigenAngle}}=nothing, mesh=nothing,
params=LMPParams())Compute electric field E, amplitude, and phase at rx.
Precomputed waveguide modes can optionally be provided as a Vector{EigenAngle}. By default modes are found with findmodes.
If mesh = nothing, use defaultmesh to generate mesh for the mode finding algorithm. This is ignored if modes is not nothing.
propagate(waveguide::SegmentedWaveguide, tx::Emitter, rx::AbstractSampler;
mesh=nothing, params=LMPParams())Compute electric field E, amplitude, and phase at rx through a SegmentedWaveguide.
If mesh = nothing, use defaultmesh to generate mesh for the mode finding algorithm.
propagate(file::AbstractString, outfile=missing; incrementalwrite=false, append=false,
mesh=nothing)Run the model scenario described by file and save the results as outfile.
If outfile = missing, the output file name will be $(file)_output.json.
LongwaveModePropagator.LMPParams — TypeLMPParams{T,T2,H <: AbstractRange{Float64}}Parameters for the LongwaveModePropagator module with defaults:
topheight::Float64 = 110e3: starting height for integration of the ionosphere reflection coefficient.earthradius::Float64 = 6369e3: Earth radius in meters.earthcurvature::Bool = true: toggle inclusion of Earth curvature in calculations. This is not supported by all functions.curvatureheight::Float64 = 50e3: reference height for Earth curvature in meters. At this height, the index of refraction is 1, and is therefore the reference height for eigenangles.approxsusceptibility::Bool = false: use a cubic interpolating spline representation ofsusceptibilityduring the integration ofdRdz.susceptibilitysplinestep::Float64 = 10.0: altitude step in meters used to build the spline representation ofsusceptibilityifapproxsusceptibility == true.grpfparams::GRPFParams = GRPFParams(100000, 1e-5, true): parameters for theGRPFcomplex root-finding algorithm.integrationparams::IntegrationParams{T} = IntegrationParams(solver=Vern7(), tolerance=1e-5): parameters passed toDifferentialEquations.jlfor integration of the ionosphere reflection coefficient.wavefieldheights::H = range(topheight, 0, length=2049): heights in meters at which wavefields will be integrated.wavefieldintegrationparams::IntegrationParams{T2} = IntegrationParams(solver=Tsit5(), tolerance=1e-6): parameters passed toDifferentialEquations.jlfor integration of the wavefields used in mode conversion. The solver cannot be lazy.radiationresistancecorrection::Bool = false: Perform a radiation resistance correction to the calculated field amplitude for transmitter antennas with altitude > 0.
The struct is created using Parameters.jl @with_kw and supports that package's instantiation capabilities, e.g.:
p = LMPParams()
p2 = LMPParams(earth_radius=6370e3)
p3 = LMPParams(p2; grpf_params=GRPFParams(100000, 1e-6, true))See also: IntegrationParams
Mode Finder
LongwaveModePropagator.findmodes — Functionfindmodes(modeequation::ModeEquation, mesh=nothing; params=LMPParams())Find EigenAngles associated with modeequation.waveguide within the domain of mesh.
mesh should be an array of complex numbers that make up the original grid over which the GRPF algorithm searches for roots of modeequation. If mesh === nothing, it is computed with defaultmesh.
There is a check for redundant modes that requires modes to be separated by at least 1 orders of magnitude greater than grpfparams.tolerance in real and/or imaginary component. For example, if grpfparams.tolerance = 1e-5, then either the real or imaginary component of each mode must be separated by at least 1e-4 from every other mode.
LongwaveModePropagator.PhysicalModeEquation — TypePhysicalModeEquation{W<:HomogeneousWaveguide} <: ModeEquationParameters for solving the physical mode equation $\det(Rg*R - I)$.
Fields:
- ea::EigenAngle
- frequency::Frequency
- waveguide::WLongwaveModePropagator.setea — Functionsetea(ea, modeequation)Return modeequation with eigenangle ea.
ea will be converted to an EigenAngle if necessary.
LongwaveModePropagator.IntegrationParams — TypeIntegrationParams{T}Parameters passed to OrdinaryDiffEq.jl during the integration of the ionosphere reflection coefficient matrix in modefinder.jl.
Fields
tolerance::Float64 = 1e-5: integrationatolandrtol.solver::T = Vern7(): aDifferentialEquations.jlsolver.dt::Float64 = 1.0: height step in meters (many methods use a variable step size).force_dtmin::Bool = false: if true, continue integration when solver reaches minimum step size.maxiters::Int = 100_000: maximum number of iterations before stopping.
EigenAngle
LongwaveModePropagator.EigenAngle — TypeEigenAnglePlane-wave propagation direction in angle θ from the vertical in radians.
Common trigonometric function values of this angle are calculated to increase performance.
Fields
θ::ComplexF64cosθ::ComplexF64sinθ::ComplexF64secθ::ComplexF64cos²θ::ComplexF64sin²θ::ComplexF64
Real θ will be automatically converted to Complex.
Technically this is an angle of incidence from the vertical, and not necessarily an eigenangle unless it is found to be associated with a propagated mode in the waveguide.
LongwaveModePropagator.attenuation — Functionattenuation(ea, frequency::Frequency; params=LMPParams())Compute attenuation of eigenangle ea at the ground.
ea will be converted to an EigenAngle if needed.
This function internally references ea to the ground.
LongwaveModePropagator.phasevelocity — Functionphasevelocity(ea)Compute the relative phase velocity $v/c$ associated with the eigenangle θ.
ea will be converted to an EigenAngle if needed.
This function internally references ea to the ground.
LongwaveModePropagator.referencetoground — Functionreferencetoground(ea::EigenAngle; params=LMPParams())Reference ea from params.curvatureheight to ground ($z = 0$).
referencetoground(sinθ; params=LMPParams())Reference $sin(θ)$ from params.curvatureheight to the ground and return $sin(θ₀)$ at the ground.
Geophysics
LongwaveModePropagator.BField — TypeBFieldBackground magnetic field vector of strength B in Tesla with direction cosines dcl, dcm, and dcn corresponding to $x$, $y$, and $z$ directions parallel, perpendicular, and up to the waveguide.
LongwaveModePropagator.Species — TypeSpeciesIonosphere constituent Species.
Fields
charge::Float64: signed species charged in Coulombs.mass::Float64: species mass in kilograms.numberdensity: a callable that returns number density in number per cubic meter as a function of height in meters.collisionfrequency: a callable that returns the collision frequency in collisions per second as a function of height in meters.
numberdensity and collisionfrequency will be converted to FunctionerWrapper types that tell the compiler that these functions will always return values of type Float64. A limited test is run to check if this is true, but otherwise it is up to the user to ensure these functions return only values of type Float64.
LongwaveModePropagator.Ground — TypeGround(ϵᵣ, σ)Isotropic earth ground characterized by a relative permittivity ϵᵣ and a conductivity σ in Siemens per meter.
LongwaveModePropagator.GROUND — ConstantDefault ground conductivity indices from LWPC.
LongwaveModePropagator.waitprofile — Functionwaitprofile(z, h′, β; cutoff_low=0, threshold=1e12)Compute the electron number density in electrons per cubic meter at altitude z in meters using Wait's exponential profile [[Wait1964]; [Thomson1993]] with parameters h′ in kilometers and β in inverse kilometers.
The profile is:
\[Nₑ = 1.43 × 10¹³ \exp(-0.15 h') \exp[(β - 0.15)(z/1000 - h')]\]
Optional arguments:
cutoff_low=0: whenzis belowcutoff_low, return zero.threshold=1e12: when density is greater thanthreshold, returnthreshold.
See also: electroncollisionfrequency, ioncollisionfrequency
References
[Wait1964]: J. R. Wait and K. P. Spies, “Characteristics of the earth-ionosphere waveguide for VLF radio waves,” U.S. National Bureau of Standards, Boulder, CO, Technical Note 300, Dec. 1964.
[Thomson1993]: N. R. Thomson, “Experimental daytime VLF ionospheric parameters,” Journal of Atmospheric and Terrestrial Physics, vol. 55, no. 2, pp. 173–184, Feb. 1993, doi: 10.1016/0021-9169(93)90122-F.
LongwaveModePropagator.electroncollisionfrequency — Functionelectroncollisionfrequency(z)Compute the electron-neutral collision frequency in collisions per second at height z in meters based on Wait's conductivity profile [[Wait1964]; [Thomson1993]].
The profile is:
\[νₑ(z) = 1.816 × 10¹¹ \exp(-0.15⋅z/1000)\]
See also: waitprofile, ioncollisionfrequency
References
[Wait1964]: J. R. Wait and K. P. Spies, “Characteristics of the earth-ionosphere waveguide for VLF radio waves,” U.S. National Bureau of Standards, Boulder, CO, Technical Note 300, Dec. 1964.
[Thomson1993]: N. R. Thomson, “Experimental daytime VLF ionospheric parameters,” Journal of Atmospheric and Terrestrial Physics, vol. 55, no. 2, pp. 173–184, Feb. 1993, doi: 10.1016/0021-9169(93)90122-F.
LongwaveModePropagator.ioncollisionfrequency — Functionioncollisionfrequency(z)Compute the ion-neutral collision frequency in collisions per second at height z in meters from [Morfitt1976].
The profile is:
\[νᵢ(z) = 4.54 × 10⁹ \exp(-0.15⋅z/1000)\]
See also: waitprofile, electroncollisionfrequency
References
[Morfitt1976]: D. G. Morfitt and C. H. Shellman, “‘MODESRCH’, an improved computer program for obtaining ELF/VLF/LF mode constants in an Earth-ionosphere waveguide,” Naval Electronics Laboratory Center, San Diego, CA, NELC/IR-77T, Oct. 1976.
Samplers
LongwaveModePropagator.Sampler — TypeSampler{T} <: AbstractSampler{T}Sampler types sample (measure) the electromagnetic field in the waveguide.
Fields
distance::T: ground distance from the transmitter in meters.fieldcomponent::Fields.Field: field component measured by theSampler.altitude::Float64: height above the ground in meters.
LongwaveModePropagator.GroundSampler — TypeGroundSampler{T} <: AbstractSampler{T}GroundSamplers are Sampler types with an altitude of zero.
Fields
distance::T: ground distance from the transmitter in meters.fieldcomponent::Fields.Field: field component measured by theGroundSampler.
LongwaveModePropagator.Fields — ModuleFieldsThis baremodule allows scoped enum-like access to electric field components Ex, Ey, and Ez.
Examples
julia> Fields.Ex
Ex::Field = 0
julia> Fields.Ey
Ey::Field = 1LongwaveModePropagator.Receiver — TypeReceiver{<:Antenna}Represent a physical receiver.
A default Receiver{VerticalDipole} is returned with Receiver().
Fields
name::String = "": receiver name.latitude::Float64 = 0.0: geographic latitude in degrees.longitude::Float64 = 0.0: geographic longitude in degrees.altitude::Float64 = 0.0: receiver altitude in meters above the ground.antenna::Antenna = VerticalDipole(): receiver antenna.
Emitters
LongwaveModePropagator.Frequency — TypeFrequencyElectromagnetic wave defined by frequency f, but also carrying angular wave frequency ω, wavenumber k, and wavelength λ, all in SI units.
LongwaveModePropagator.Transmitter — TypeTransmitter{A<:Antenna} <: EmitterTypical ground-based Transmitter.
Fields
name::String: transmitter name.latitude::Float64: transmitter geographic latitude in degrees.longitude::Float64: transmitter geographic longitude in degrees.antenna::Antenna: transmitter antenna.frequency::Frequency: transmit frequency.power::Float64: transmit power in Watts.
LongwaveModePropagator.Dipole — TypeDipoleDipole antenna with arbitrary orientation described by azimuth_angle and inclination_angle from vertical.
See also: AbstractDipole
LongwaveModePropagator.VerticalDipole — TypeVerticalDipoleDipole antenna with inclination angle $γ = 0$ from the vertical.
LongwaveModePropagator.HorizontalDipole — TypeHorizontalDipoleDipole antenna with inclination angle $γ = π/2$ from the vertical and azimuth_angle orientation $ϕ$ where:
| ϕ [rad] | orientation |
|---|---|
| 0 | end fire |
| π/2 | broadside |
LongwaveModePropagator.inclination — Functioninclination(d::AbstractDipole)Return the inclination (dip) angle $γ$ of d measured from vertical.
See also: AbstractDipole
LongwaveModePropagator.azimuth — Functionazimuth(d::AbstractDipole)Return the azimuth angle $ϕ$ of d measured towards positive y direction from x.
See also: AbstractDipole
Waveguides
LongwaveModePropagator.HomogeneousWaveguide — TypeHomogeneousWaveguide{S} <: WaveguideDefines a homogeneous segment of waveguide.
Fields
bfield::BField: background magnetic field.species::S: ionosphere constituents.ground::Ground: waveguide ground.distance::Float64: distance from theEmitterat the start of the segment in meters.
LongwaveModePropagator.SegmentedWaveguide — TypeSegmentedWaveguide{T<:Vector{<:Waveguide}} <: WaveguideA collection of Waveguides make up an inhomogeneous segmented waveguide.
IO
LongwaveModePropagator.ExponentialInput — TypeExponentialInputType for Wait and Spies (1964) exponential ionosphere profiles defined by waitprofile. The electroncollisionfrequency is used for the electron-neutral collision frequency profile.
- The electron density profile begins at 40 km altitude and extends to 110 km.
- The transmitter power is 1 kW.
Fields
name::Stringdescription::Stringdatetime::DateTimesegment_ranges::Vector{Float64}: distance from transmitter to the beginning of eachHomogeneousWaveguidesegment in meters.hprimes::Vector{Float64}: Wait's $h′$ parameter for eachHomogeneousWaveguidesegment.betas::Vector{Float64}: Wait's $β$ parameter for eachHomogeneousWaveguidesegment.b_mags::Vector{Float64}: magnetic field magnitude for eachHomogeneousWaveguidesegment.b_dips::Vector{Float64}: magnetic field dip angles in radians for eachHomogeneousWaveguidesegment.b_azs::Vector{Float64}: magnetic field azimuth in radians "east" of the propagation direction for eachHomogeneousWaveguidesegment.ground_sigmas::Vector{Float64}: ground conductivity in Siemens per meter for eachHomogeneousWaveguidesegment.ground_epsrs::Vector{Int}: ground relative permittivity for eachHomogeneousWaveguidesegment.frequency::Float64: transmitter frequency in Hertz.output_ranges::Vector{Float64}: distances from the transmitter at which the field will be calculated.
LongwaveModePropagator.TableInput — TypeTableInput <: InputFields
name::Stringdescription::Stringdatetime::DateTimesegment_ranges::Vector{Float64}: distance from transmitter to the beginning of eachHomogeneousWaveguidesegment in meters.altitude::Vector{Float64}: altitude above ground in meters for which thedensityandcollision_frequencyprofiles are specified.density::Vector{Float64}: electron density at eachaltitudein $m⁻³$.collision_frequency::Vector{Float64}: electron-ion collision frequency at eachaltitudein $s⁻¹$.b_dips::Vector{Float64}: magnetic field dip angles in radians for eachHomogeneousWaveguidesegment.b_azs::Vector{Float64}: magnetic field azimuth in radians "east" of the propagation direction for eachHomogeneousWaveguidesegment.ground_sigmas::Vector{Float64}: ground conductivity in Siemens per meter for eachHomogeneousWaveguidesegment.ground_epsrs::Vector{Int}: ground relative permittivity for eachHomogeneousWaveguidesegment.frequency::Float64: transmitter frequency in Hertz.output_ranges::Vector{Float64}: distances from the transmitter at which the field will be calculated.
LongwaveModePropagator.BatchInput — TypeBatchInput{T} <: InputA collection of inputs with a batch name, description, and datetime.
LongwaveModePropagator.BasicOutput — TypeBasicOutput <: OutputFields
name::Stringdescription::Stringdatetime::DateTimeoutput_ranges::Vector{Float64}amplitude::Vector{Float64}phase::Vector{Float64}
LongwaveModePropagator.BatchOutput — TypeBatchOutput{T} <: OutputA collection of outputs with a batch name, description, and datetime.
See also: BatchInput