Internals
Documentation of all const
s, struct
s, function
s, etc.
LongwaveModePropagator.AbstractDipole
— TypeAbstractDipole <: Antenna
Abstract type for dipole antennas.
Subtypes of AbstractDipole
have an orientation azimuth_angle
($ϕ$) and inclination_angle
($γ$) where
[rad] | ϕ | γ |
---|---|---|
0 | end fire | vertical |
π/2 | broadside | horizontal |
Angles describe the orientation of the antenna, not the radiation pattern.
LongwaveModePropagator.Antenna
— TypeAntenna
Abstract supertype for all antenna-like types.
LongwaveModePropagator.bookerquartic
— Functionbookerquartic(ea::EigenAngle, M)
bookerquartic(T::TMatrix)
Compute roots q
and the coefficients B
of the Booker quartic described by the susceptibility tensor M
or T
matrix.
References
[Budden1988]: K. G. Budden, “The propagation of radio waves: the theory of radio waves of low power in the ionosphere and magnetosphere,” First paperback edition. New York: Cambridge University Press, 1988.
LongwaveModePropagator.bookerreflection
— Functionbookerreflection(ea::EigenAngle, M::SMatrix{3,3})
bookerreflection(ea::EigenAngle, e)
Compute the ionosphere reflection coefficient matrix for a sharply bounded ionosphere from 4×2 wavefields matrix e
or the susceptibility matrix M
.
The ionosphere reflection coefficient matrix is computed from a ratio of the downgoing to upgoing plane waves in the free space beneath the ionosphere [Budden1988] pg. 307. These are obtained from the two upgoing characteristic waves found from the Booker quartic. Each make up a column of e
.
\[R = \begin{pmatrix} Ce₁[4] - e₁[1] & Ce₂[4] - e₂[1] \\ -Ce₁[2] + e₁[3] & -Ce₂[2] + e₂[3] \end{pmatrix} \begin{pmatrix} Ce₁[4] + e₁[1] & Ce₂[4] + e₂[1] \\ -Ce₁[2] - e₁[3] & -Ce₂[2] - e₂[3] \end{pmatrix}^{-1}\]
The reflection coefficient matrix for the sharply bounded case is commonly used as a starting solution for integration of the reflection coefficient matrix through the ionosphere.
References
[Budden1988]: K. G. Budden, “The propagation of radio waves: the theory of radio waves of low power in the ionosphere and magnetosphere,” First paperback edition. New York: Cambridge University Press, 1988.
Extended help
The set of horizontal field components $e = (Ex, -Ey, Z₀Hx, Z₀Hy)ᵀ$ can be separated into an upgoing and downgoing wave, each of which is generally elliptically polarized. A ratio of the amplitudes of these two waves give a reflection coefficient, except it would only apply for an incident wave of that particular elliptical polarization. However, the first set of fields can be linearly combined with a second independent solution for the fields, which will generally have a different elliptical polarization than the first. Two linear combinations of the two sets of fields are formed with unit amplitude, linearly polarized incident waves. The reflected waves then give the components $R₁₁$, $R₂₁$ or $R₁₂$, $R₂₂$ for the incident wave in the plane of incidence and perpendicular to it, respectively [Budden1988] pg 552.
The process for determining the reflection coefficient requires resolving the two sets of fields $e₁$ and $e₂$ into the four linearly polarized vacuum modes. The layer of vacuum can be assumed to be so thin that it does not affect the fields. There will be two upgoing waves and two downgoing waves, each which has one $E$ and one $H$ in the plane of incidence. If $f₁, f₂, f₃, f₄$ are the complex amplitudes of the four component waves, then in matrix notation $e = Lf$ where $L$ is the appropriate transformation matrix.
For $e₁$ and $e₂$, we can find the corresponding vectors $f1$ and $f2$ by $f1 = L⁻¹e₁$, $f2 = L⁻¹e₂$ where the two column vectors are partitioned such that $f1 = (u1, d1)ᵀ$ and $f2 = (u2, d2)ᵀ$ for upgoing and downgoing 2-element vectors $u$ and $d$. From the definition of the reflection coefficient $R$, $d = Ru$. Letting $U = (u1, u2)$, $D = (d1, d2)$, then $D = RU$ and the reflection coefficient is $R = DU⁻¹$. Because the reflection coefficient matrix is a ratio of fields, either $e₁$ and/or $e₂$ can be independently multiplied by an arbitrary constant and the value of $R$ is unaffected.
This function directly computes $D$ and $U$ and solves for $R$ using the right division operator R = D/U
.
For additional details, see [Budden1988], chapter 18, section 7.
LongwaveModePropagator.bookerreflection
— Methodbookerreflection(ea::EigenAngle, M, ::Dθ)
Compute the ionosphere reflection coefficient matrix $R$ for a sharply bounded ionosphere with susceptibility tensor M
, as well as its derivative $dR/dθ$ returned as the tuple (R, dR)
.
LongwaveModePropagator.bookerwavefields
— Functionbookerwavefields(ea::EigenAngle, M)
bookerwavefields(T::TMatrix)
Compute the two-column wavefields matrix e
from the ionosphere with susceptibility tensor M
or T
matrix for the two upgoing wavefields.
The first column of e
is the evanescent wave and the second is the travelling wave.
\[e = \begin{pmatrix} Ex₁ & Ex₂ \\ -Ey₁ & -Ey₂ \\ Hx₁ & Hx₂ \\ Hy₁ & Hy₂ \end{pmatrix}\]
This function solves the eigenvalue problem $Te = qe$. First, the Booker quartic is solved for the roots q
. Then they are sorted so that the roots associated with the two upgoing waves can be selected. e
is solved as the eigenvectors for the two q
s. An analytical solution is used where e[2,:] = 1
.
LongwaveModePropagator.bookerwavefields
— Methodbookerwavefields(ea::EigenAngle, M, ::Dθ)
bookerwavefields(T::TMatrix, dT, ::Dθ)
Compute the two-column wavefields matrix e
as well as its derivative with respect to $θ$, returned as a tuple (e, de)
for the ionosphere with susceptibility tensor M
or T
matrix and its derivative with respect to $θ$, dT
.
LongwaveModePropagator.dbookerquartic
— Functiondbookerquartic(ea::EigenAngle, M, q, B)
dbookerquartic(T::TMatrix, dT, q, B)
Compute derivative dq
of the Booker quartic roots q
with respect to $θ$ for the ionosphere described by susceptibility tensor M
or T
matrix.
LongwaveModePropagator.sortquarticroots!
— Methodsortquarticroots!(q)
Sort array of quartic roots q
in place such that the first two correspond to upgoing waves and the latter two correspond to downgoing waves.
General locations of the four roots of the Booker quartic on the complex plane corresponding to the:
1) upgoing evanescent wave
2) upgoing travelling wave
3) downgoing evanescent wave
4) downgoing travelling wave
Im
3 . |
|
|
. 4 |
------------------Re
| . 2
|
|
| . 1
Based on [Pitteway1965] fig. 5.
References
[Pitteway1965]: M. L. V. Pitteway, “The numerical calculation of wave-fields, reflexion coefficients and polarizations for long radio waves in the lower ionosphere. I.,” Phil. Trans. R. Soc. Lond. A, vol. 257, no. 1079, pp. 219–241, Mar. 1965, doi: 10.1098/rsta.1965.0004.
LongwaveModePropagator.upgoing
— Methodupgoing(q)
Calculate the absolute angle of q
in radians from 315°×π/180 on the complex plane. Smaller values indicate upgoing waves.
Base.isless
— Methodisless(x::EigenAngle, y::EigenAngle)
Calls Base.isless
for complex EigenAngle
by real
and then imag
component.
By defining isless
, calling sort
on EigenAngle
s sorts by real component first and imaginary component second.
LongwaveModePropagator.isdetached
— Methodisdetached(ea::EigenAngle, frequency::Frequency; params=LMPParams())
Return true
if this is likely an earth detached (whispering gallery) mode according to the criteria in [Pappert1981] eq. 1 with the additional criteria that the frequency be above 100 kHz.
References
[Pappert1981]: R. A. Pappert, “LF daytime earth ionosphere waveguide calculations,” Naval Ocean Systems Center, San Diego, CA, NOSC/TR-647, Jan. 1981. [Online]. Available: https://apps.dtic.mil/docs/citations/ADA096098.
LongwaveModePropagator.AirborneTransmitter
— TypeAirborneTransmitter{A<:Antenna} <: Emitter
Similar to Transmitter
except it also contains an altitude
field in meters.
LongwaveModePropagator.Emitter
— TypeEmitter
Abstract type for types that energize the waveguide with electromagnetic energy.
LongwaveModePropagator.dip
— Methoddip(b::BField)
Return the dip angle in radians from b
.
LongwaveModePropagator.gyrofrequency
— Methodgyrofrequency(q, m, B)
gyrofrequency(s::Species, b::BField)
Return the signed gyrofrequency $Ω = qB/m$.
LongwaveModePropagator.magnetoionicparameters
— Methodmagnetoionicparameters(z, frequency::Frequency, bfield::BField, species::Species)
Compute the magnetoionic parameters X
, Y
, and Z
for height z
.
\[X = N e² / (ϵ₀ m ω²) Y = e B / (m ω) Z = ν / ω\]
References
[Budden1955a]: K. G. Budden, “The numerical solution of differential equations governing reflexion of long radio waves from the ionosphere,” Proc. R. Soc. Lond. A, vol. 227, no. 1171, pp. 516–537, Feb. 1955, pp. 517.
[Budden1988]: K. G. Budden, “The propagation of radio waves: the theory of radio waves of low power in the ionosphere and magnetosphere,” First paperback edition. New York: Cambridge University Press, 1988, pp. 39.
[Ratcliffe1959]: J. A. Ratcliffe, "The magneto-ionic theory & its applications to the ionosphere," Cambridge University Press, 1959.
LongwaveModePropagator.plasmafrequency
— Methodplasmafrequency(n, q, m)
plasmafrequency(z, s::Species)
Return the plasma frequency $ωₚ = √(nq²/(ϵ₀m))$ for a ``cold'' plasma.
LongwaveModePropagator.waitsparameter
— Methodwaitsparameter(z, frequency::Frequency, bfield::BField, species::Species)
Compute Wait's conductivity parameter Wr
.
\[ωᵣ = Xω²/ν = ωₚ²/ν\]
LongwaveModePropagator.Input
— TypeInput
Abstract supertype for structs carrying information to be input to the model.
LongwaveModePropagator.Output
— TypeOutput
Abstract supertype for structs containing information to be output from the model.
LongwaveModePropagator.buildrun
— Functionbuildrun(s::ExponentialInput; mesh=nothing, unwrap=true, params=LMPParams())
buildrun(s::TableInput; mesh=nothing, unwrap=true, params=LMPParams())
buildrun(s::BatchInput; mesh=nothing, unwrap=true, params=LMPParams())
Build LMP structs from an Input
and run LMP
.
For TableInput
s, a FritschButland monotonic interpolation is performed over density
and collision_frequency
.
LongwaveModePropagator.buildrunsave
— Methodbuildrunsave(outfile, s::BatchInput; append=false, mesh=nothing, unwrap=true, params=LMPParams())
Similar to buildrun
, except it saves results into outfile
as s
is processed.
If append=true
, this function parses outfile
for preexisting results and only runs the remaining scenarios in s
. Otherwise, a new BatchOutput
is created.
LongwaveModePropagator.buildwaveguide
— Methodbuildwaveguide(s::ExponentialInput, i)
Return HomogeneousWaveguide
from the i
th entry in each field of s
.
LongwaveModePropagator.buildwaveguide
— Methodbuildwaveguide(s::TableInput, i)
Return HomogeneousWaveguide
from the i
th entry in each field of s
with a FritschButland monotonic interpolation over density
and collision_frequency
.
Outside of s.altitude
the nearest s.density
or s.collision_frequency
is used.
LongwaveModePropagator.iscomplete
— Methodiscomplete(s)
Return true
if input or output struct s
is completely defined, otherwise return false
.
LongwaveModePropagator.parse
— Methodparse(file)
Parse a JSON file compatible with Input
or Output
types.
LongwaveModePropagator.validlengths
— Functionvalidlengths(s)
Check if field lengths of input s
match their number of segments.
LongwaveModePropagator.BOTTOMHEIGHT
— ConstantBottom of reflection coefficient integration. LongwaveModePropagator.jl
assumes BOTTOMHEIGHT = 0.0
.
LongwaveModePropagator.Dθ
— TypeIndicates derivative with respect to eigenangle $θ$.
LongwaveModePropagator.ModeEquation
— TypeModeEquation
Functions can dispatch on subtypes of ModeEquation
, although currently only PhysicalModeEquation
is supported.
Future work might include the ModifiedModeEquation
of [Morfitt1976].
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. [Online]. Available: http://www.dtic.mil/docs/citations/ADA032573.
LongwaveModePropagator.susceptibility
— Methodsusceptibility(altitude, frequency, bfield, species; params=LMPParams())
susceptibility(altitude, frequency, w::HomogeneousWaveguide; params=LMPParams())
susceptibility(altitude, me::ModeEquation; params=LMPParams())
Compute the ionosphere susceptibility tensor M
as a SMatrix{3,3}
using species.numberdensity
and species.collisionfrequency
at altitude
.
Multiple species can be passed as an iterable. Use a tuple
of Species
, rather than a Vector
, for better performance.
If params.earthcurvature == true
, M
includes a first order correction for earth curvature by means of a fictitious refractive index [Pappert1967].
The susceptibility matrix is calculated from the constitutive relations presented in [Ratcliffe1959]. This includes the effect of earth's magnetic field vector and collisional damping on electron motion.
The tensor is:
\[M = -\frac{X}{U(U²-Y²)} \begin{pmatrix} U² - x²Y² & -izUY - xyY² & iyUY - xzY² \\ izUY - xyY² & U² - y²Y² & -ixUY - yzY² \\ -iyUY - xzY² & ixUY - yzY² & U² - z²Y² \end{pmatrix}\]
where $X = ωₚ²/ω²$, $Y = |ωₕ/ω|$, $Z = ν/ω$, and $U = 1 - iZ$. The earth curvature correction subtracts $2/Rₑ*(H - altitude)$ from the diagonal of $M$ where $H$ is params.curvatureheight
.
References
[Pappert1967]: R. A. Pappert, E. E. Gossard, and I. J. Rothmuller, “A numerical investigation of classical approximations used in VLF propagation,” Radio Science, vol. 2, no. 4, pp. 387–400, Apr. 1967, doi: 10.1002/rds196724387.
[Ratcliffe1959]: J. A. Ratcliffe, "The magneto-ionic theory & its applications to the ionosphere," Cambridge University Press, 1959.
LongwaveModePropagator.susceptibilityspline
— Methodsusceptibilityspline(frequency, bfield, species; params=LMPParams())
susceptibilityspline(frequency, w::HomogeneousWaveguide; params=LMPParams())
susceptibilityspline(me::ModeEquation; params=LMPParams())
Construct a cubic interpolating spline of susceptibility
and return the callable Interpolations
type.
params.susceptibilitysplinestep
is the altitude step in meters used to construct the spline between BOTTOMHEIGHT
and params.topheight
.
LongwaveModePropagator.modeconversion
— Methodmodeconversion(previous_wavefields::Wavefields, wavefields::Wavefields,
adjoint_wavefields::Wavefields; params=LMPParams())
Compute the mode conversion matrix a
from the modes associated with previous_wavefields
to modes associated with wavefields
and its adjoint_wavefields
.
This is used in the approach known as full mode conversion [Pappert1972b].
References
[Pappert1972b]: R. A. Pappert and R. R. Smith, “Orthogonality of VLF height gains in the earth ionosphere waveguide,” Radio Science, vol. 7, no. 2, pp. 275–278, 1972, doi: 10.1029/RS007i002p00275.
LongwaveModePropagator.dRdz
— FunctiondRdz(R, modeequation, z, susceptibilityfcn=z->susceptibility(z, modeequation; params=LMPParams()))
Compute the differential of the reflection matrix R
, $dR/dz$, at height z
. susceptibilityfcn
is a function returning the ionosphere susceptibility at height z
.
Following the Budden formalism for the reflection of an (obliquely) incident plane wave from a horizontally stratified ionosphere [Budden1955a], the differential of the reflection matrix R
with height z
can be described by
\[dR/dz = k/(2i)⋅(W₂₁ + W₂₂R - RW₁₁ - RW₁₂R)\]
Integrating $dR/dz$ downwards through the ionosphere gives the reflection matrix $R$ for the ionosphere as if it were a sharp boundary at the stopping level with free space below.
References
[Budden1955a]: K. G. Budden, “The numerical solution of differential equations governing reflexion of long radio waves from the ionosphere,” Proc. R. Soc. Lond. A, vol. 227, no. 1171, pp. 516–537, Feb. 1955.
LongwaveModePropagator.dRdθdz
— MethoddRdθdz(RdRdθ, p, z)
Compute the differential $dR/dθ/dz$ at height z
returned as an SMatrix{4,2}
with $dR/dz$ in the first 2 rows and $dR/dθ/dz$ in the bottom 2 rows.
p
is a tuple containing instances (PhysicalModeEquation(), LMPParams())
.
LongwaveModePropagator.defaultmesh
— Methoddefaultmesh(frequency; rmax=deg2rad(89.9), imax=deg2rad(0.0),
Δr_coarse=deg2rad(0.5), Δr_fine=deg2rad(0.1),
rtransition=deg2rad(75.0), itransition=deg2rad(-1.5),
meshshape="auto")
Generate vector of complex coordinates (radians) to be used by GRPF in the search for waveguide modes.
rmin
is the lower bound of the real axis and imin
is the lower bound of the imaginary axis.
The value of frequency
sets different default behavior:
frequency | meshshape = "auto" | rmin | imin | resolution |
---|---|---|---|---|
< 12 kHz | "rectanglemesh" | deg2rad(1) | deg2rad(-89) | Δr_coarse |
≥ 12 kHz | "trianglemesh" | deg2rad(30) | deg2rad(-10) | variable |
At frequencies at or above 12 kHz the mesh spacing in the upper right corner of the domain with real values above rtransition
and imaginary values above itransition
is Δr_fine
and is Δr_coarse
everywhere else.
If meshshape = "rectanglemesh"
, the full lower right quadrant of the complex plane is searched for modes. If meshshape = "trianglemesh"
, the lower right diagonal of the lower right quadrant of the complex plane is excluded from the mesh. Eigenangles at VLF (~12 kHz) and higher frequencies are not typically in the lower right diagonal.
See also: findmodes
LongwaveModePropagator.dmodalequation
— Methoddmodalequation(R, dR, Rg, dRg)
Compute the derivative of the determinantal mode equation with respect to $θ$.
LongwaveModePropagator.dwmatrix
— Methoddwmatrix(ea::EigenAngle, T, dT)
Compute the four submatrix elements of $dW/dθ$ returned as the tuple (dW₁₁, dW₂₁, dW₁₂, dW₂₂)
from the ionosphere with T
matrix and its derivative with respect to $θ$, dT
.
LongwaveModePropagator.fresnelreflection
— Functionfresnelreflection(ea::EigenAngle, ground::Ground, frequency::Frequency)
fresnelreflection(m::PhysicalModeEquation)
Compute the Fresnel reflection coefficient matrix for the ground-freespace interface at the ground.
LongwaveModePropagator.fresnelreflection
— Methodfresnelreflection(ea::EigenAngle, ground::Ground, frequency::Frequency, ::Dθ)
fresnelreflection(m::PhysicalModeEquation, ::Dθ)
Compute the Fresnel reflection coefficient matrix for the ground as well as its derivative with respect to $θ$ returned as the tuple (Rg, dRg)
.
LongwaveModePropagator.integratedreflection
— Methodintegratedreflection(modeequation::PhysicalModeEquation, ::Dθ; params=LMPParams())
Compute $R$ and $dR/dθ$ as an SMatrix{4,2}
with $R$ in rows (1, 2) and $dR/dθ$ in rows (3, 4).
The params.integrationparams.tolerance
is hardcoded to 1e-10
in this version of the function.
LongwaveModePropagator.integratedreflection
— Methodintegratedreflection(modeequation::PhysicalModeEquation;
params=LMPParams(), susceptibilityfcn=z->susceptibility(z, modeequation; params=params))
Integrate $dR/dz$ downward through the ionosphere described by modeequation
from params.topheight
, returning the ionosphere reflection coefficient R
at the ground. susceptibilityfcn
is a function returning the ionosphere susceptibility tensor as a function of altitude z
in meters.
params.integrationparams
are passed to DifferentialEquations.jl
.
LongwaveModePropagator.modalequation
— Methodmodalequation(R, Rg)
Compute the determinental mode equation $det(Rg R - I)$ given reflection coefficients R
and Rg
.
A propagating waveguide mode requires that a wave, having reflected from the ionosphere and then the ground, must be identical with the original upgoing wave. This criteria is met at roots of the mode equation [Budden1962].
References
[Budden1962]: K. G. Budden and N. F. Mott, “The influence of the earth’s magnetic field on radio propagation by wave-guide modes,” Proceedings of the Royal Society of London. Series A. Mathematical and Physical Sciences, vol. 265, no. 1323, pp. 538–553, Feb. 1962.
LongwaveModePropagator.solvedmodalequation
— Methodsolvedmodalequation(θ, modeequation::PhysicalModeEquation; params=LMPParams())
Set θ
for modeequation
and then solve the derivative of the mode equation with respect to θ
.
LongwaveModePropagator.solvedmodalequation
— Methodsolvedmodalequation(modeequation::PhysicalModeEquation; params=LMPParams())
Compute the derivative of the modal equation with respect to $θ$ returned as the tuple (dF, R, Rg)
for the ionosphere and ground reflection coefficients.
LongwaveModePropagator.solvemodalequation
— Methodsolvemodalequation(θ, modeequation::PhysicalModeEquation;
params=LMPParams(), susceptibilityfcn=z->susceptibility(z, modeequation; params=params))
Set θ
for modeequation
and then solve the modal equation.
LongwaveModePropagator.solvemodalequation
— Methodsolvemodalequation(modeequation::PhysicalModeEquation;
params=LMPParams(), susceptibilityfcn=z->susceptibility(z, modeequation; params=params))
Compute the ionosphere and ground reflection coefficients and return the value of the determinental modal equation associated with modeequation
. susceptibilityfcn
is a function that returns the ionosphere susceptibility as a function of altitude z
in meters.
See also: solvedmodalequation
LongwaveModePropagator.trianglemesh
— Methodtrianglemesh(zbl, ztr, Δr)
Generate initial mesh node coordinates for a right triangle domain from complex coordinate zbl
in the bottom left and ztr
in the top right with initial mesh step Δr
.
zbl
and ztr
are located on the complex plane at the locations marked in the diagram below:
im
|
-re ---|------ ztr
| /
| /
| /
| /
zbl
LongwaveModePropagator.wmatrix
— Methodwmatrix(ea::EigenAngle, T)
Compute the four submatrix elements of W
used in the equation $dR/dz$ from the ionosphere with T
matrix returned as a tuple (W₁₁, W₂₁, W₁₂, W₂₂)
.
Following Budden's formalism for the reflection matrix of a plane wave obliquely incident on the ionosphere [Budden1955a], the wave below the ionosphere can be resolved into upgoing and downgoing waves of elliptical polarization, each of whose components are themselves resolved into a component with the electric field in the plane of propagation and a component perpendicular to the plane of propagation. The total field can be written in matrix form as $e = Lf$ where $L$ is a 4×4 matrix that simply selects and specifies the incident angle of the components and $f$ is a column matrix of the complex amplitudes of the component waves. By inversion, $f = L⁻¹e$ and its derivative with respect to height $z$ is $f′ = -iL⁻¹TLf = -½iWf$. Then $W = 2L⁻¹TL$ describes the change in amplitude of the upgoing and downgoing component waves.
$W$ is also known as $S$ in many texts.
References
[Budden1955a]: K. G. Budden, “The numerical solution of differential equations governing reflexion of long radio waves from the ionosphere,” Proc. R. Soc. Lond. A, vol. 227, no. 1171, pp. 516–537, Feb. 1955.
LongwaveModePropagator.ExcitationFactor
— TypeExcitationFactor{T,T2}
Constants used in calculating excitation factors and height gains.
Fields
F₁::T
: height gain constant. See [Pappert1976].F₂::T
F₃::T
F₄::T
h₁0::T
: first modified Hankel function of order 1/3 at the ground.h₂0::T
: second modified Hankel function of order 1/3 at the ground.EyHy::T
: polarization ratio $Ey/Hy$, derived from reflection coefficients (or $T$s).Rg::T2
: ground reflection coefficient matrix.
References
[Pappert1976]: R. A. Pappert and L. R. Shockey, “Simplified VLF/LF mode conversion program with allowance for elevated, arbitrarily oriented electric dipole antennas,” Naval Electronics Laboratory Center, San Diego, CA, Interim Report 771, Oct. 1976. [Online]. Available: http://archive.org/details/DTIC_ADA033412.
LongwaveModePropagator.excitationfactor
— Methodexcitationfactor(ea, dFdθ, R, Rg, efconstants::ExcitationFactor; params=LMPParams())
Compute excitation factors for the $Hy$ field at the emitter returned as the tuple (λv, λb, λe)
for vertical, broadside, and end-on dipoles. dFdθ
is the derivative of the modal equation with respect to $θ$.
The excitation factor describes how efficiently the field component can be excited in the waveguide.
This function most closely follows the approach taken in [Pappert1983], which makes use of $T$ (different from TMatrix
) rather than $τ$. From the total $Hy$ excitation factor (the sum product of the λ
s with the antenna orientation terms), the excitation factor for electric fields can be found as:
- $λz = -S₀λ$
- $λx = EyHy⋅λ$
- $λy = -λ$
This function assumes that reflection coefficients R
and Rg
are referenced to $d = z = 0$.
References
[Morfitt1980]: D. G. Morfitt, “‘Simplified’ VLF/LF mode conversion computer programs: GRNDMC and ARBNMC,” Naval Ocean Systems Center, San Diego, CA, NOSC/TR-514, Jan. 1980. [Online]. Available: http://www.dtic.mil/docs/citations/ADA082695.
[Pappert1983]: R. A. Pappert, L. R. Hitney, and J. A. Ferguson, “ELF/VLF (Extremely Low Frequency/Very Low Frequency) long path pulse program for antennas of arbitrary elevation and orientation,” Naval Ocean Systems Center, San Diego, CA, NOSC/TR-891, Aug. 1983. [Online]. Available: http://www.dtic.mil/docs/citations/ADA133876.
[Pappert1986]: R. A. Pappert and J. A. Ferguson, “VLF/LF mode conversion model calculations for air to air transmissions in the earth-ionosphere waveguide,” Radio Sci., vol. 21, no. 4, pp. 551–558, Jul. 1986, doi: 10.1029/RS021i004p00551.
LongwaveModePropagator.excitationfactorconstants
— Methodexcitationfactorconstants(ea₀, R, Rg, frequency, ground; params=LMPParams())
Return an ExcitationFactor
struct used in calculating height-gain functions and excitation factors where eigenangle ea₀
is referenced to the ground.
This function assumes that reflection coefficients R
and Rg
are referenced to $d = z = 0$.
References
[Pappert1976]: R. A. Pappert and L. R. Shockey, “Simplified VLF/LF mode conversion program with allowance for elevated, arbitrarily oriented electric dipole antennas,” Naval Electronics Laboratory Center, San Diego, CA, Interim Report 771, Oct. 1976. [Online]. Available: http://archive.org/details/DTIC_ADA033412.
[Ferguson1980]: J. A. Ferguson and F. P. Snyder, “Approximate VLF/LF waveguide mode conversion model: Computer applications: FASTMC and BUMP,” Naval Ocean Systems Center, San Diego, CA, NOSC-TD-400, Nov. 1980. [Online]. Available: http://www.dtic.mil/docs/citations/ADA096240.
[Morfitt1980]: D. G. Morfitt, “‘Simplified’ VLF/LF mode conversion computer programs: GRNDMC and ARBNMC,” Naval Ocean Systems Center, San Diego, CA, NOSC/TR-514, Jan. 1980. [Online]. Available: http://www.dtic.mil/docs/citations/ADA082695.
LongwaveModePropagator.fieldsum
— Methodfieldsum(modes, waveguide::HomogeneousWaveguide, tx::Emitter, rx::AbstractSampler;
params=LMPParams())
Compute the complex electric field by summing modes
in waveguide
for emitter tx
at sampler rx
.
This function always returns all three electric field components, regardless of the value of rx.fieldcomponent
.
References
[Morfitt1980]: D. G. Morfitt, “‘Simplified’ VLF/LF mode conversion computer programs: GRNDMC and ARBNMC,” Naval Ocean Systems Center, San Diego, CA, NOSC/TR-514, Jan. 1980. [Online]. Available: http://www.dtic.mil/docs/citations/ADA082695.
[Pappert1983]: R. A. Pappert, L. R. Hitney, and J. A. Ferguson, “ELF/VLF (Extremely Low Frequency/Very Low Frequency) long path pulse program for antennas of arbitrary elevation and orientation,” Naval Ocean Systems Center, San Diego, CA, NOSC/TR-891, Aug. 1983. [Online]. Available: http://www.dtic.mil/docs/citations/ADA133876.
LongwaveModePropagator.fieldsum
— Methodfieldsum(waveguide::SegmentedWaveguide, wavefields_vec, adjwavefields_vec, tx::Emitter,
rx::AbstractSampler; params=LMPParams())
LongwaveModePropagator.heightgains
— Methodheightgains(z, ea₀, frequency, efconstants::ExcitationFactor; params=LMPParams())
Compute height gain functions at height z
returned as the tuple (fz, fy, fx)
where eigenangle ea₀
is referenced to the ground.
fz
is the height gain for the vertical electric field component $Ez$.fy
is the height gain for the transverse electric field component $Ey$.fx
is the height gain for the horizontal electric field component $Ex$.
[Pappert1983]
This function assumes that reflection coefficients are referenced to $d = z = 0$.
See also: excitationfactorconstants
References
[Pappert1983]: R. A. Pappert, L. R. Hitney, and J. A. Ferguson, “ELF/VLF (Extremely Low Frequency/Very Low Frequency) long path pulse program for antennas of arbitrary elevation and orientation,” Naval Ocean Systems Center, San Diego, CA, NOSC/TR-891, Aug. 1983. [Online]. Available: http://www.dtic.mil/docs/citations/ADA133876.
[Pappert1986]: R. A. Pappert and J. A. Ferguson, “VLF/LF mode conversion model calculations for air to air transmissions in the earth-ionosphere waveguide,” Radio Sci., vol. 21, no. 4, pp. 551–558, Jul. 1986, doi: 10.1029/RS021i004p00551.
LongwaveModePropagator.modeterms
— Methodmodeterms(modeequation, tx::Emitter, rx::AbstractSampler; params=LMPParams())
Compute tx
and rx
height-gain and excitation factor products and ExcitationFactor
constants returned as the tuple (txterm, rxterm)
.
The returned txterm
is:
\[λ_v \cos(γ) f_z(zₜ) + λ_b \sin(γ)\sin(ϕ) f_y(zₜ) + λ_e \sin(γ)\cos(ϕ) f_z(zₜ)\]
and rxterm
is the height-gain function $f(zᵣ)$ appropriate for rx.fieldcomponent
:
fieldcomponent | $f(zᵣ)$ |
---|---|
$z$ | $-S₀⋅f_z$ |
$y$ | $EyHy⋅f_y$ |
$x$ | $-f_x$ |
References
[Pappert1976]: R. A. Pappert and L. R. Shockey, “Simplified VLF/LF mode conversion program with allowance for elevated, arbitrarily oriented electric dipole antennas,” Naval Electronics Laboratory Center, San Diego, CA, Interim Report 771, Oct. 1976. [Online]. Available: http://archive.org/details/DTIC_ADA033412.
[Morfitt1980]: D. G. Morfitt, “‘Simplified’ VLF/LF mode conversion computer programs: GRNDMC and ARBNMC,” Naval Ocean Systems Center, San Diego, CA, NOSC/TR-514, Jan. 1980. [Online]. Available: http://www.dtic.mil/docs/citations/ADA082695.
LongwaveModePropagator.radiationresistance
— Methodradiationresistance(k, Cγ, zt)
Calculate radiation resistance correction for transmitting antenna elevated above the ground. Based on [Pappert1986] below, is derived from the time-averaged Poynting vector, which uses total E and H fields calculated assuming a point dipole over perfectly reflecting ground.
If a point dipole at height z is radiating known power Pz, the power that should be input to LMP is $P/Pz = 2/f(kz,γ)$. Also note that E ∝ √P, hence the square root below.
References
[Pappert1986]: R. A. Pappert, “Radiation resistance of thin antennas of arbitrary elevation and configuration over perfectly conducting ground.,” Naval Ocean Systems Center, San Diego, CA, Technical Report 1112, Jun. 1986. Accessed: Mar. 10, 2024. [Online]. Available: https://apps.dtic.mil/sti/citations/ADA170945
LongwaveModePropagator.AbstractSampler
— TypeAbstractSampler
Abstract supertype for sampling fields in the waveguide.
Subtypes of AbstractSampler have a position in the guide and a field component.
LongwaveModePropagator.TMatrix
— TypeTMatrix <: SMatrix{4, 4, T}
A custom SMatrix
subtype that represents the T
matrix from [Clemmow1954] that is semi-sparse. The form of the matrix is:
┌ ┐
| T₁₁ T₁₂ 0 T₁₄ |
| 0 0 1 0 |
| T₃₁ T₃₂ 0 T₃₄ |
| T₄₁ T₄₂ 0 T₄₄ |
└ ┘
TMatrix
implements efficient matrix-vector multiplication and other math based on its special form.
See also: tmatrix
References
[Clemmow1954]: P. C. Clemmow and J. Heading, “Coupled forms of the differential equations governing radio propagation in the ionosphere,” Mathematical Proceedings of the Cambridge Philosophical Society, vol. 50, no. 2, pp. 319–333, Apr. 1954.
LongwaveModePropagator.dtmatrix
— Methoddtmatrix(ea::EigenAngle, M)
Compute a dense SMatrix
with the derivative of T
with respect to θ
.
See also: tmatrix
LongwaveModePropagator.tmatrix
— Methodtmatrix(ea::EigenAngle, M)
Compute the matrix T
as a TMatrix
for the differential equations of a wave propagating at angle ea
in an ionosphere with susceptibility M
.
Clemmow and Heading derived the T
matrix from Maxwell's equations for an electromagnetic wave in the anisotropic, collisional cold plasma of the ionosphere in a coordinate frame where $z$ is upward, propagation is directed obliquely in the $x$-$z$ plane and invariance is assumed in $y$. For the four characteristic wave components $e = (Ex, -Ey, Z₀Hx, Z₀Hy)ᵀ$, the differential equations are $de/dz = -ikTe$.
See also: susceptibility
, dtmatrix
References
[Budden1955a]: K. G. Budden, “The numerical solution of differential equations governing reflexion of long radio waves from the ionosphere,” Proc. R. Soc. Lond. A, vol. 227, no. 1171, pp. 516–537, Feb. 1955.
[Clemmow1954]: P. C. Clemmow and J. Heading, “Coupled forms of the differential equations governing radio propagation in the ionosphere,” Mathematical Proceedings of the Cambridge Philosophical Society, vol. 50, no. 2, pp. 319–333, Apr. 1954.
LongwaveModePropagator.amplitude
— Methodamplitude(e)
Compute field amplitude in dB.
LongwaveModePropagator.amplitudephase
— Methodamplitudephase(e)
Compute field amplitude in dB and phase in radians and return as (amplitude
, phase
).
LongwaveModePropagator.pow23
— Functionpow23(x)
Efficiently compute $x^(2/3)$.
LongwaveModePropagator.unwrap!
— Functionunwrap!(x)
Unwrap a phase vector x
in radians in-place.
LongwaveModePropagator.ScaleRecord
— TypeScaleRecord
Struct for saving wavefield scaling information used in a callback during Pitteway integration of wavefields [Pitteway1965].
Fields
z::Float64
: current height in meters.e::SMatrix{4,2,Complex{Float64},8}
: wavefield matrix.ortho_scalar::Complex{Float64}
: scalar for Gram-Schmidt orthogonalization.e1_scalar::Float64
: scalar for wavefield vector 1.e2_scalar::Float64
: scalar for wavefield vector 2.
LongwaveModePropagator.WavefieldIntegrationParams
— TypeWavefieldIntegrationParams{S,T,T2,H}
Parameters passed to Pitteway integration of wavefields [Pitteway1965].
Fields
z::Float64
: height in meters.bottomz::Float64
: bottom height of integration in meters.ortho_scalar::Complex{Float64}
: scalar for Gram-Schmidt orthogonalization.e1_scalar::Float64
: scalar for wavefield vector 1.e2_scalar::Float64
: scalar for wavefield vector 2.ea::EigenAngle
: wavefield eigenangle.frequency::Frequency
: electromagentic wave frequency in Hz.bfield::BField
: Earth's magnetic field in Tesla.species::S
: species in the ionosphere.params::LMPParams{T,T2,H}
: module-wide parameters.
LongwaveModePropagator.WavefieldIntegrationParams
— MethodWavefieldIntegrationParams(topheight, ea, frequency, bfield, species, params)
Initialize a WavefieldIntegrationParams
for downward Pitteway scaled integration [Pitteway1965].
Automatically set values are:
bottomz = BOTTOMHEIGHT
ortho_scalar = zero(ComplexF64)
e1_scalar = one(Float64)
e2_scalar = one(Float64)
LongwaveModePropagator.Wavefields
— TypeWavefields{H<:AbstractRange}
Struct containing wavefield components for each mode of eas
at heights
.
The six electromagnetic field components are stored in the v
field of the struct and accessed as a Matrix
of SVector{6,ComplexF64}
corresponding to [height, mode]
.
LongwaveModePropagator.Wavefields
— MethodWavefields(heights, eas::Vector{EigenAngle})
A constructor for initializing Wavefields
with an appropriately sized undef
Matrix given eigenangles eas
and heights
.
Base.:==
— Method==(A::Wavefields, B::Wavefields)
Check for equality ==
between A.v
and B.v
, heights
, and eas
.
LongwaveModePropagator.boundaryscalars
— Functionboundaryscalars(R, Rg, e1, e2, isotropic::Bool=false)
boundaryscalars(R, Rg, e, isotropic::Bool=false)
Compute coefficients (b1, b2)
required to sum the two wavefields vectors e1
and e2
or both columns of e
for the total wavefield at the ground as $e = b1*e1 + b2*e2$.
The b1
and b2
that satisfy the waveguide boundary conditions are only valid for true eigenangles of the waveguide.
This function assumes that the reflection coefficients R
and Rg
and the wavefield vectors e1
, e2
are at the ground.
These boundary conditions
LongwaveModePropagator.calculate_wavefields!
— Methodcalculate_wavefields!(wavefields, adjoint_wavefields, frequency, waveguide,
adjoint_waveguide; params=LMPParams())
Compute fields of wavefields
in-place scaled to satisfy the waveguide
boundary conditions.
This function implements the method of integrating wavefields suggested by [Pitteway1965].
References
[Pitteway1965]: M. L. V. Pitteway, “The numerical calculation of wave-fields, reflexion coefficients and polarizations for long radio waves in the lower ionosphere. I.,” Phil. Trans. R. Soc. Lond. A, vol. 257, no. 1079, pp. 219–241, Mar. 1965, doi: 10.1098/rsta.1965.0004.
LongwaveModePropagator.dedz
— Methoddedz(e, p, z)
Compute derivative of field components vector e
at height z
.
The parameter p
should be a WavefieldIntegrationParams
.
LongwaveModePropagator.dedz
— Methoddedz(e, k, T::Tmatrix)
Compute $de/dz = -ikTe$ where $e = (Ex, -Ey, Z₀Hx, Z₀Hy)ᵀ$.
LongwaveModePropagator.fieldstrengths!
— Methodfieldstrengths!(EH, zs, me::ModeEquation; params=LMPParams())
fieldstrengths!(EH, zs, ea::EigenAngle, frequency::Frequency, bfield::BField,
species, ground::Ground; params=LMPParams())
Compute $(Ex, Ey, Ez, Hx, Hy, Hz)ᵀ$ wavefields vectors as elements of EH
by fullwave integration at each height in zs
.
The wavefields are scaled to satisfy the waveguide boundary conditions, which is only valid at solutions of the mode equation.
LongwaveModePropagator.fieldstrengths
— Methodfieldstrengths(zs, me::ModeEquation; params=LMPParams())
Preallocate vector of wavefields EH
, then call fieldstrengths!
and return EH
.
Each element of EH
is an SVector
of $ex, ey, ez, hx, hy, hz$.
LongwaveModePropagator.integratewavefields
— Methodintegratewavefields(zs, ea::EigenAngle, frequency::Frequency, bfield::BField,
species; params=LMPParams(), unscale=true)
Compute wavefields vectors e
at zs
by downward integration over heights zs
.
params.wavefieldintegrationparams
is used by this function rather than params.integrationparams
.
unscale
corrects scaling of wavefields with unscalewavefields
.
LongwaveModePropagator.savevalues
— Methodsavevalues(u, t, integrator)
Return a ScaleRecord
from u
, t
, and integrator
.
LongwaveModePropagator.scale!
— Methodscale!(integrator)
Apply wavefield scaling with scalewavefields
to the integrator.
LongwaveModePropagator.scalewavefields
— Functionscalewavefields(e1, e2)
scalewavefields(e::SMatrix{4,2})
Orthonormalize the vectors e1
and e2
or the columns of e
, and also return the scaling terms a
, e1_scale_val
, and e2_scale_val
applied to the original vectors.
This first applies Gram-Schmidt orthogonalization and then scales the vectors so they each have length 1, i.e. norm(e1) == norm(e2) == 1
. This is the technique suggested by [Pitteway1965] to counter numerical swamping during integration of wavefields.
References
[Pitteway1965]: M. L. V. Pitteway, “The numerical calculation of wave-fields, reflexion coefficients and polarizations for long radio waves in the lower ionosphere. I.,” Phil. Trans. R. Soc. Lond. A, vol. 257, no. 1079, pp. 219–241, Mar. 1965, doi: 10.1098/rsta.1965.0004.
LongwaveModePropagator.scalingcondition
— Methodscalingcondition(e, z, int)
Return true
if wavefields should be scaled, otherwise false
for wavefields e
at height z
and integrator int
.
Wavefields should be scaled if any component of real(e)
or imag(e)
are >= 1
. In addition, force scaling at z <= bottomz
to ensure initial upgoing wave is unit amplitude.
LongwaveModePropagator.unscalewavefields!
— Methodunscalewavefields!(e, saved_values::SavedValues)
Unscale the integrated wavefields, a vector of fields at each integration step e
, in place.
See also: unscalewavefields
LongwaveModePropagator.unscalewavefields
— Methodunscalewavefields(saved_values::SavedValues)
Return the unscaled integrated wavefields originally scaled by scalewavefields
.
The bottom level does not get unscaled. We reference the higher levels to the bottom. The level above the bottom level is additionally scaled by the amount that was applied to originally get from this level down to the bottom level. The next level up (2 above the bottom level) is scaled by the amount applied to the next level and then the bottom level, i.e. we keep track of a cumulative correction on the way back up.
Assumes fields have been scaled by scalewavefields
during integration.
See also: unscalewavefields!
LongwaveModePropagator.Waveguide
— TypeWaveguide
A waveguide propagation path with a background BField
, ionosphere Species
, and Ground
.
LongwaveModePropagator.adjoint
— Methodadjoint(w::HomogeneousWaveguide)
Return w
with an adjoint BField
having an x
component of opposite sign.