API Reference
Publicly Exported Functions
Base.complex — Methodcomplex(t::Exi)Return a vector of complex excitation amplitudes.
Base.isapprox — Methodisapprox(c1::Cut, c2::Cut; kwargs...) -> tf::BoolCompare two Cut objects for approximate equality.
Compares most fields for perfect equality except text and evec. The text fields are not compared at all, and the evec fields (Matrix types) are compared for approximate equality using isapprox.
Base.maximum — Methodmaximum(cut::Cut) -> maxEReturn the maximum amplitude for any polarization component stored in the Cut object.
Base.:+ — Method+(s1::Surface, s2::Surface)Return a new Surface whose $z$ values are the sum of those of s1 and s2. s1 and s2 must have identical x and y fields.
Base.:- — Method-(s1::Surface, s2::Surface)Return a new Surface whose $z$ values are the difference of those of s1 and s2. s1 and s2 must have identical x and y fields.
TicraUtilities.amplitude_db — Functionamplitude_db(c::Cut, ipol::Int)
amplitude_db(c::Cut, polsymb::Symbol = :copol)
amplitude_db(c::Cut, polstr::String = "copol")Return a matrix of amplitudes in dB for some choice of polarization in the cut. Legal values for ipol are 1, 2, or 3, the latter only being legal if there are three polarization components present in the cut. Legal values for polstr are "copol" (the default) and "xpol". Capitalization is not significant. Legal values for polsymb are :copol and :xpol. Again, capitalization is not significant. Copol is defined as the polarization with maximum amplitude at θ = ϕ = 0.
TicraUtilities.amplitude_db — Methodamplitude_db(t::Exi)Return a vector of excitation amplitudes in dB
TicraUtilities.asym2sym — Methodasym2sym(cut::Cut) -> cut2Convert an asymmetrical (in θ) Cut to a new symmetrical Cut. An asymmetrical cut begins at θ = 0, while a symmetrical cut covers equal extents of negative and positive angles. If the input cut is already symmetrical, then return a copy of this input as the output.
TicraUtilities.convert_cut! — Methodconvert_cut!(cut::Cut, icomp::Integer)Convert cut to a new polarization basis determined by icomp. Legal values for icomp and their meanings:
- 1 => Eθ and Eϕ
- 2 => ERHCP and ELHCP
- 3 => Eh and Ev (Ludwig 3 co and cx)
TicraUtilities.convert_cut — Methodcut2 = convert_cut(cut::Cut, icomp::Integer)Make a copy of cut and convert it to a new polarization basis determined by icomp. Legal values for icomp and their meanings:
- 1 => Eθ and Eϕ
- 2 => ERHCP and ELHCP
- 3 => Eh and Ev (Ludwig 3 co and cx)
TicraUtilities.cut2sph — Methodcut2sph(cut::Cut; keywords...) -> sph::SPHQPartition
cut2sph(cuts::AbstractVector{Cut}; keywords...) -> sphs::Vector{SPHQPartition}
cut2sph(cutfile::AbstractString; kwargs...) -> s::SPHQPartitionConvert a Cut object to a SPHQPartition using recursive FFT/IFFT methods from the Hansen 1988 book "Spherical Near-Field Antenna Measurements.
The single positional input argument can be either a string containing the name of a Ticra-compatible, spherical polar cut file, or the returned value of type Cut (or a vector of Cut objects) that results from reading such a file with read_cutfile. The output of this function can be passed to write_sphfile to create a Ticra-compatible file of Q-type spherical wave coefficients. If the input cut contains 3 polarization slots, the third slot will be discarded before performing the transformation.
If the input cuts extend in θ only to θ₀ < 180°, then it will be assumed that the fields are identically zero for θ₀ < θ ≤ 180°.
Keyword Arguments (and their default values)
mmax=1000: An upper limit for them(azimuthal) mode index to be included. The actual limit will be set tomin(Nϕ÷2, mmax)for oddNϕ, andmin(Nϕ÷2-1, mmax)for evenNϕ, whereNϕis the number of ϕ = constant polar cuts in the cut object.nmax=1000: An upper limit for then(polar) mode index to be included. The actual limit will be the lesser ofnmaxandNθ-1whereNθis the number of θ values included in each ϕ = constant polar cut.pwrtol=0.0: The power tolerance. Spherical modes are included until the excluded modes' power is less thanpwrtoltimes the total modal power. A zero or negative value precludes removal of any modes.
TicraUtilities.discard_pol3 — Methoddiscard_pol3(cut::Cut) -> cut2::CutCreate a copy of cut but discard the third polarization slot, if present.
TicraUtilities.eh2bor1cut — Methodeh2bor1cut(theta, fe, fh; kwargs...) -> cut::CutCreate a Cut object for a "BOR1" horn from its E-plane and H-plane patterns.
A "BOR₁" horn is circularly symmetric and contains only TE₁ₙ and TM₁ₙ waveguide modes in its radiating aperture. It's radiated far field can therefore be expressed in terms of the E-plane and H-plane patterns it radiates when excited for linear polarization.
Positional Input Arguments
theta: A vector or range (anAbstractVector) of θ values (in degrees) at which the cut pattern should be evaluated. The first element ofthetamust be 0, and the entries must be equally spaced, as in arangeobject.fe,fh: The E-plane and H-plane patterns, resp. These are either bothAbstractVectors of the same length astheta, or both functions which take a single input θ (in degrees) and return the respective patterns evaluated at that angle.
Keyword Arguments
pol: Defines the manner in which the horn is assumed to be excited, and the polarization basis selected for use in the outputCut.polis aStringorSymboltaking one of the values (capitalization is not significant):- "l3v" or
:l3v: (the default value) The horn is excited for "vertical" (y-directed) linear polarization and the far field is expressed as Ludwig-3 components. - "l3h" or
:l3h: The horn is excited for "horizontal" (x-directed) linear polarization and the far field is expressed as Ludwig-3 components. - "rhcp" or
:rhcp: The horn is excited for RHCP (right-hand circular polarization) and the far field is expressed as RHCP and LHCP components. - "lhcp" or
:lhcp: The horn is excited for LHCP (left-hand circular polarization) and the far field is expressed as RHCP and LHCP components.
Cutobject will contain eight (four) cuts, spaced every 45° (90°).- "l3v" or
xpd: The crosspol level in dB < 0. Defaults to-Inf(negative infinity). If finite, then in addition to the specified polarization, a crosspolarized contribution will be added to the cut, as if the horn is fed by an imperfect feed network with the specified crosspol level.xpphase: The phase (in degrees) of the crosspol contribution whose amplitude is specified byxpd.
TicraUtilities.get_ampdb — Methodget_ampdb(t::Exi)Return a vector of excitation amplitudes in dB
TicraUtilities.get_att — Methodget_att(s::Station)Return att, the vector of attenuation values wrt nadir in dB.
TicraUtilities.get_class — Methodget_class(tep::TEPperiodic)Return the class name string for the TEP object.
TicraUtilities.get_evec — Methodget_evec(c::Cut, ipol::Integer)Return the ntheta × nphi matrix of complex numbers stored in polarization slot ipol of the cut. ipol must be positive and less than or equal to get_ncomp(cut).
TicraUtilities.get_evec — Methodget_evec(c::Cut)Return the ntheta × nphi matrix of complex field vectors stored in the cut. Each element of the returned matrix will be either a 2-vector or 3-vector, depending on the number of field components stored in the cut.
TicraUtilities.get_freqs — Methodget_freqs(tep::TEPperiodic)Return the vector of frequencies. The elements of the vector will be Unitful quantities, each of which may have different frequency units.
TicraUtilities.get_goal — Methodget_goal(s::Station)Return goal, the vector of optimization goals in dB.
TicraUtilities.get_header — Methodget_header(t::Exi)Return a vector of excitation header strings.
TicraUtilities.get_icomp — Methodget_icomp(c::Cut)Return icomp, polarization parameter. 1 for Eθ and Eφ, 2 for RHCP and LHCP, 3 for h and v (Ludwig 3).
TicraUtilities.get_icut — Methodget_icut(c::Cut)Return icut, the control parameter of the cut. 1 for a polar cut, 2 for a conical cut. TicraUtilities currently accommodates only icut == 1, wherein each cut is for a constant ϕ value.
TicraUtilities.get_id — Methodget_id(s::Station)Return id, the vector of station ID strings.
TicraUtilities.get_ids — Methodget_ids(t::Exi)Return a vector of excitation ID strings.
TicraUtilities.get_idstrg — Methodget_idstrg(s::SPHQPartition)Return the idstrg string of s.
TicraUtilities.get_ipol — Methodget_ipol(s::Station)Return ipol, the vector of station polarization codes.
TicraUtilities.get_mmax — Methodget_mmax(s::SPHQPartition)Return the integer mmax (maximum m index) associated with s.
TicraUtilities.get_name — Methodget_name(tep::TEPperiodic)Return the object name string for the TEP object.
TicraUtilities.get_name — Methodget_name(obj::TorObj) -> name::StringReturn the object name.
TicraUtilities.get_ncomp — Methodget_ncomp(c::Cut)Return ncomp, the number of polarization components stored in the cut.
TicraUtilities.get_nmax — Methodget_nmax(s::SPHQPartition)Return the integer nmax (maximum n index) associated with s.
TicraUtilities.get_nphi — Methodget_nphi(s::SPHQPartition)Return the integer nphi (number of phi points) associated with s.
TicraUtilities.get_npoint — Methodget_npoint(s::Station)Return npoint, the number of stations.
TicraUtilities.get_nthe — Methodget_nthe(s::SPHQPartition)Return the integer nthe (number of theta points) associated with s.
TicraUtilities.get_objtype — Methodget_objtype(obj::TorObj) -> objtype::StringReturn the Ticra object type.
TicraUtilities.get_phi — Methodget_phi(c::Cut)Return the phi values (degrees) stored in the cut. The returned value will be some type of AbstractRange object.
TicraUtilities.get_phi — Methodget_phi(tep::TEP)Return the phi values (degrees) stored in the TEP object. The returned value will be some type of AbstractRange object.
TicraUtilities.get_phsdeg — Methodget_phsdeg(t::Exi)Return a vector of excitation phases in degrees.
TicraUtilities.get_powerms — Methodget_powerms(s::SPHQPartition)Return the offset vector powerms. The vector has axes (0:mmax) (i.e. it can be indexed with integers m ranging from 0 to mmax) and its mth element contains the total power (one-half the sum of the magnitude squared) of all modes with ±m as the m index.
TicraUtilities.get_prgtag — Methodget_prgtag(s::SPHQPartition)Return the prgtag string of s.
TicraUtilities.get_propname — Methodget_propname(obj::TorObj) -> propname::Vector{String}Return the vector of object property names.
TicraUtilities.get_propval — Methodget_propval(obj::TorObj) -> propval::Vector{String}Return the vector of object property values.
TicraUtilities.get_qsmns — Methodget_qsmns(s::SPHQPartition)Return the offset array qsmns::OffsetArray{ComplexF64, 3} that contains the Q coefficients, assuming exp(jωt) time variation and using the Ticra normalization convention. The array has axes (1:2, -mmax:mmax, 1:nmax), meaning that some of the entries (those with n < abs(m)) are always zero.
TicraUtilities.get_rot — Methodget_rot(s::Station)Return rot, the vector of station polarization rotation angles in degrees.
TicraUtilities.get_sff — Methodget_sff(tep::TEP)Return the reflection array for front incidence. The array has size (2, 2, length(theta), length(phi) for a TEPscatter object and (2, 2, length(theta), length(phi), length(freqs)) for a TEPperiodic object.
TicraUtilities.get_sfr — Methodget_sfr(tep::TEP)Return the transmission array for rear incidence. The array has size (2, 2, length(theta), length(phi) for a TEPscatter object and (2, 2, length(theta), length(phi), length(freqs)) for a TEPperiodic object.
TicraUtilities.get_srf — Methodget_srf(tep::TEP)Return the transmission array for front incidence. The array has size (2, 2, length(theta), length(phi) for a TEPscatter object and (2, 2, length(theta), length(phi), length(freqs)) for a TEPperiodic object.
TicraUtilities.get_srr — Methodget_srr(tep::TEP)Return the reflection array for rear incidence. The array has size (2, 2, length(theta), length(phi) for a TEPscatter object and (2, 2, length(theta), length(phi), length(freqs)) for a TEPperiodic object.
TicraUtilities.get_t4 — Methodget_t4(s::SPHQPartition)Return the string t4 associated with s.
TicraUtilities.get_t5 — Methodget_t5(s::SPHQPartition)Return the string t5 associated with s.
TicraUtilities.get_t6 — Methodget_t6(s::SPHQPartition)Return the string t6 associated with s.
TicraUtilities.get_t7 — Methodget_t7(s::SPHQPartition)Return the string t7 associated with s.
TicraUtilities.get_t8 — Methodget_t8(s::SPHQPartition)Return the string t8 associated with s.
TicraUtilities.get_text — Methodget_text(c::Cut)Return a vector of strings containing the cut identification text.
TicraUtilities.get_text — Methodget_text(c::Surface) -> t::StringReturn the text string associated with surface c.
TicraUtilities.get_theta — Methodget_theta(c::Cut)Return the theta values (in degrees) stored in the cut. The returned value will be some type of AbstractRange object.
TicraUtilities.get_theta — Methodget_theta(tep::TEP)Return the theta values (degrees) stored in the TEP object. The returned value will be some type of AbstractRange object.
TicraUtilities.get_title — Methodget_title(tep::TEPscatter)Return the title string for the TEP object.
TicraUtilities.get_u — Methodget_u(s::Station)Return u, the vector of unitless station direction cosines along $x$.
TicraUtilities.get_v — Methodget_v(s::Station)Return v, the vector of unitless station direction cosines along $y$.
TicraUtilities.get_weight — Methodget_weight(s::Station)Return weight, the vector of station optimization weights.
TicraUtilities.get_x — Methodget_x(c::Surface) -> x::AbstractRangeReturn the range of x values for surface c.
TicraUtilities.get_y — Methodget_y(c::Surface) -> y::AbstractRangeReturn the range of y values for surface c.
TicraUtilities.get_z — Methodget_z(c::Surface) -> z::Matrix{Float64}Return the matrix of z values for surface c. z[i,j] is the $z$ value sampled at $(x, y)$ coordinates $($x[i], y[j]$)$.
TicraUtilities.issym — Methodissym(cut::Cut) -> tf::BoolCheck if a cut has a symmetrical range in theta that includes zero.
TicraUtilities.maximum_db — Methodmaximum_db(cut::Cut)Return the maximum amplitude in dB for any polarization component stored in the Cut object.
TicraUtilities.normalize! — Functionnormalize!(cut::Cut, totpower=4π)Normalize a Cut object so it's total power is totpower (which defaults to 4π). The default value results in field magnitude squared being numerically equal to directivity.
TicraUtilities.parse_tor_struct — Methodparse_tor_struct(s::AbstractString) -> t::NamedTupleThe input string is an entry in a TorObj propval field, consisting of, e.g., "struct(status: automatic, x: 1 in, y:23.4 mm)". This function parses out the Ticra struct and returns a NamedTuple t whose propertynames are the Ticra struct field names. The values of t are either strings or Unitful quantities. In the example given, t.status == "automatic", t.x == 1.0u"inch", and t.y == 23.4u"mm".
TicraUtilities.phase_deg — Functionphase_deg(c::Cut, ipol::Int)
phase_deg(c::Cut, polsymb::Symbol)
phase_deg(c::Cut, polstr::String = "copol")Return a matrix of phases in degrees for some choice of polarization in the cut. Legal values for ipol are 1, 2, or 3. Legal values for polstr are "copol" (the default) and "xpol" (capitalization is not significant). Legal values of polsymb are :copol and :xpol. Again, capitalization is not significant. The returned value is a Matrix{Float64} with length(get_theta(cut)) rows and length(get_phi(cut)) columns.
TicraUtilities.phase_deg — Methodphase_deg(t::Exi)Return a vector of excitation phases in degrees.
TicraUtilities.phscen — Function(x,y,z0,z90) = phscen(cut, fghz=11.80285268; min_dropoff=-10)Estimate phase center for a cut file or Cut object using NSI least-squares algorithm.
The four outputs are estimates of the location of the phase center relative to the origin used in recording the data in the cut. If fghz is passed in as an argument, the values will be expressed in units of inches. Otherwise the lengths will be normalized to wavelengths. Note that z0 and z90 are the phi = 0ᵒ and phi = 90ᵒ plane estimates of the phase center location. In determining the phase center locations, only field values with magnitudes in dB greater than min_dropoff relative to the peak field are considered.
TicraUtilities.power — Functionpower(cut::Cut, θmax=180)Compute the total radiated power in a Cut object. If only a single phi value is included in the cut, then assume no azimuthal variation. The integration in the θ direction will be computed over the limits from 0 to min(θmax, last(get_theta(cut))). θmax should be specified in degrees.
TicraUtilities.read_cutfile — Methodread_cutfile(fname) -> Vector{Cut}Read data from a possibly multi-frequency Ticra-compatible cut file.
Return a single Cut struct or a vector Cut structs. Each element of the returned vector corresponds to a particular operating frequency partition in the file. If there is only a single partition in the file, then instead of returning a 1-element vector, the single object of type Cut is returned as a scalar.
TicraUtilities.read_exifile — Methodread_exifile(filename::AbstractString)Reads contents of a Ticra-compatible excitation file and returns a Exi object.
TicraUtilities.read_sphfile — Methodread_sphfile(fname) -> Vector{SPHQPartition}Read the SPH coefficients from a Q-type spherical wave expansion file.
In the process of reading the data, the coefficients in the file (Q′) are conjugated and then multiplied by the factor √(8π) to achieve Ticra-standard normalization. Each element of the returned vector corresponds to a particular operating frequency partition in the file. If there is only a single partition in the file, then instead of returning a 1-element vector, the single element of type SPHQPartition is returned as a scalar.
TicraUtilities.read_stationfile — Methodread_stationfile(stationfile) -> Vector{Station}Read the contents of a Ticra optimization station file. The returned value is a vector of length npart, where npart is the number of partitions in the file.
TicraUtilities.read_surface — Methodread_surface(fname::AbstractString) -> s::SurfaceRead a Ticra-compatible surface file and return a Surface object.
TicraUtilities.read_tepfile — Methodread_tepfile(filename::AbstractString)Read a TICRA-compatible "Tabulated Electrical Properties" (TEP) file. The file may be in either the original format (scattering surface) originated by GRASP8 or the newer format for periodic unit cells created by the QUPES program. The return value depends on the type of TEP file encountered:
- Old-style scattering surface: For a TEP file containing data for a single frequency, the return value will be a scalar object of type
TEPscatter. If the file contains data for multiple frequencies, the return value will be an object of typeVector{TEPscatter}, with one element for each frequency. - New-style periodic unit cell: The return value will be an object of type
TEPperiodic. Note that the file may not contain swept geometrical parameters.
TicraUtilities.read_torfile — Methodread_torfile(torfile::AbstractString)Return a vector of TorObj objects found in a TOR file
TicraUtilities.sor_efficiency — Methodsor_efficiency(cut; F, D, Oc, pol=:matched, dz=0.0) -> (;ηₗₒₛₛ, ηₛₚ, ηᵢₗ, ηₚₕ, ηₚₒₗ, ηₓ)Compute boresight directivity efficiency of a parabolic single-offset reflector using the feed pattern specified by cut.
Positional input arguments:
cut: Either a string containing the name of a Ticra-compatible, spherical, polar, asymmetric cut file, or aCutobject as returned by theread_cutfilefunction.
Keyword input arguments:
F,D, andOc: The single-offset reflector focal length, aperture diameter, and center offset, respectively. These may be expressed in any convenient length units, so long as they are consistent.pol: ASymbolhaving one of the values:L3h,:L3v,:RHCP,:LHCP, or :max (any variations in terms of capitalization are acceptable). The first two denote Ludwig 3 horizontal (x) and vertical (y) polarizations, the second two denote the two senses of circular polarization, and:max(the default) uses the polarization among the 4 previously listed that that has the maximum field amplitude. Polarization efficiency and boresight polarization mismatch of the secondary pattern will be computed relative to the polarization specified by this argument. Note that for CP, this will be orthogonal to the primary copol.dz: This is a signed distance along the zfeed direction, measured in wavelengths. It allows for repositioning the feed in an attempt to locate the feed phase center at the reflector focal point. Suppose that the feed's phase center is located 0.42 wavelengths inside the horn aperture. The horn origin should ideally be positioned closer to the reflector, so a positive valuedz = 0.42would be specified to indicate that the horn has been repositioned in this manner.
Return value: A NamedTuple with the following fields:
;ηₗₒₛₛ, ηₛₚ, ηᵢₗ, ηₚₕ, ηₚₒₗ, ηₓ
ηₗₒₛₛ: Feed loss efficiency (ratio of radiated power to 4π).ηₛₚ: Spillover efficiency (a real number between 0 and 1).ηᵢₗ: Illumination (amplitude taper) efficiency (a real number between 0 and 1).ηₚₕ: Phase error efficiency (a real number between 0 and 1).ηₚₒₗ: Polarization efficiency (a real number between 0 and 1).ηₓ: Boresight polarization mismatch loss efficiency (a real number between 0 and 1)
TicraUtilities.sph2cut — Methodsph2cut(sphfile::AbstractString; theta, phi, ipol) -> cut::Cut
sph2cut(sph:SPHQPartition; theta, phi, ipol) -> cut::Cut
sph2cut(sphs:Vector{SPHQPartition}; theta, phi, ipol) -> cuts::Vector{Cut}Convert a set of Q-type spherical wave modal coefficients to far-field electric field values, returned as a Cut object.
The single positional input argument can be either a string containing the name of a Ticra-compatible Q-type spherical wave file, or the returned value from reading such a file with read_sphfile.
Optional Keyword Arguments:
theta: An abstract range denoting the desired polar angles (colattitudes) in degrees at which the field should be evaluated. If an empty range is provided (the default), then the values will be determined automatically by examining the modal content insph.phi: An abstract range denoting the desired azimuthal angles in degrees at which the field should be evaluated. If an empty range is provided (the default), then the values will be determined automatically by examining the modal content insph.ipol: An integer in the range 0:3 denoting the desired polarization decomposition for the computed field values. The meanings of these values are:- 0: (Default value) Choose the basis among choices 2 or 3 that produces the largest peak copol magnitude.
- 1: Use a (θ̂, ϕ̂) basis.
- 2: Use a (R̂, L̂) (i.e., circular polarization) basis.
- 3: Use a (ĥ, v̂) (Ludwig 3) basis.
Usage Example
cut = sph2cut("testfile.sph"; phi=0:5:355, theta=0:1:180, ipol=2)TicraUtilities.sym2asym — Methodsym2asym(cut::Cut) -> cut2Convert a symmetrical (in θ) Cut to a new asymmetrical Cut. An asymmetrical cut begins at θ = 0, while a symmetrical cut covers equal extents of negative and positive angles. If the input cut is already asymmetrical, then return a copy of this input as the output.
TicraUtilities.symsqueeze — Methodsymsqueeze(cut::Cut) -> cut2::CutSqueeze a symmetric cut with ϕ values equally distributed in [0,360) into [0, 180).
The ϕ values must define nonredundant cuts. I.e., cut.phi must be equivalent to the vector [0, Δϕ, 2Δϕ, 3Δϕ, ..., (360 - Δϕ)], with length(cut.phi) being an odd number, so that 180 is excluded.
cut2 will have Δϕ/2 as its ϕ increment, and it's maximum ϕ value will be 180 - Δϕ/2.
TicraUtilities.tepp2s — Methodtepp2s(tep::TEPperiodic, d; title="")Convert a periodic unit cell TEP object (of type TEPperiodic) to a scattering surface TEP object (of type TEPscatter), or to a vector of such objects if tep contains more than a single frequency.
d is the is the distance between the front and rear reference planes, a Unitful length quantity. The title keyword argument is used for the title field of the output object. If it is left at its default empty value, then it will be replaced by "TEPscatter object created by tepp2s".
Usage Example
using Unitful: @u_str
d = 2.3u"mm"
tep_scatter = tepp2s(tep_periodic, d)Extended help
Input argument d is required because TEPperiodic uses phase reference planes located at the actual front and rear surfaces of the unit cell, while TEPscatter uses the front surface only as the phase reference plane for both front and rear incidence. Thus d is required to compute the necessary phase correction for rear surface incidence reflection and for both front and rear surface incidence transmission. This phase correction is in addition to sign changes needed for some of the coefficients.
TicraUtilities.teps2p — Functionteps2p(tep::TEPscatter, d, freq; name="tep_periodic", class="created_by_teps2p")Convert a scattering TEP object (of type TEPscatter) to a periodic unit cell TEP object (of type TEPperiodic).
d is the is the distance between the front and rear reference planes, a Unitful quantity. freq is the frequency, also a Unitful quantity. tep and freq may be both scalars or both vectors of the same length. In the latter case, each entry corresponds to a specific frequency.
Usage Example
using Unitful: @u_str
freqs = [1.2u"GHz", 2u"GHz"] # Assumes teps is a vector of 2 TEPscatter objects
d = 2.3u"mm"
tep_periodic = teps2p(teps, d, freqs)Extended help
Input argument freq is required because the frequencies are part of the data stored in a TEPperiodic object, but not in a TEPscatter object. Additionally, both freq and d arguments are required because TEPperiodic uses phase reference planes located at the actual front and rear surfaces of the unit cell, while TEPscatter uses the front surface only as the phase reference plane for both front and rear incidence. Thus, these arguments are required to compute the necessary phase correction for rear surface incidence reflection and for both front and rear surface incidence transmission. This phase correction is in addition to sign changes needed for some of the coefficients.
TicraUtilities.write_cutfile — Methodwrite_cutfile(fname::AbstractString, cut::Cut, title::AbstractString="Cut file created by write_cutfile")
write_cutfile(fname::AbstractString, cuts::AbstractVector{Cut}, title::AbstractString="Cut file created by write_cutfile")Write Cut cut data to a Ticra-compatible cut file.
TicraUtilities.write_exifile — Methodwrite_exifile(filename::AbstractString, t::Exi)Create a Ticra-compatible excitation file from a Exi object.
TicraUtilities.write_sphfile — Methodwrite_sphfile(fname, qs::Vector{SPHQPartition})
write_sphfile(fname, qs::SPHQPartition)Write SPH coefficients to a Q-type spherical wave expansion file.
Prior to writing the data into the file, the input coefficients (Q) are conjugated and then multiplied by the factor 1/sqrt(8π) to become Q′ and achieve consistency with Ticra-standard normalization.
TicraUtilities.write_stationfile — Methodwrite_stationfile(stationfile::AbstractString, stdat::Station)
write_stationfile(stationfile::AbstractString, stdat::AbstractVector{Station})Write a Ticra POS4-compatible optimization station file. Here, when stdat is a vector, its elements are the partitions in the station file.
TicraUtilities.write_surface — Methodwrite_surface(fname::AbstractString, sfc::Surface)Write a Surface object to a Ticra-compatible surface file.
TicraUtilities.write_tepfile — Functionwrite_tepfile(filename::AbstractString, tep::TEP)
write_tepfile(filename::AbstractString, tep::Vector{TEPscatter})Write a TICRA-compatible "Tabulated Electrical Properties" (TEP) file. If type(tep) == TEPscatter or type(tep) == Vector{TEPscatter}, then the file will be written in the original format (scattering surface) introduced by GRASP8. If type(tep) == TEPperiodic then the file will be written in the newer format for periodic unit cells introduced by the QUPES program.
TicraUtilities.write_torfile — Methodwrite_torfile(fname::AbstractString, obj::TorObj; modestr = "w")
write_torfile(fname::AbstractString, objs::Vector{TorObj}; modestr = "w")Write one or multiple TorObj objects to a Ticra-compatible TOR (Ticra Object Repository) file. With the default value of modestr = "w", the file is created if it doesn't already exist, or, if it already exists, it is truncated to zero length before writing the TOR object(s). By setting the value of modestr to "a" one can append the object(s) to an existing file.
Publicly Exported Types
TicraUtilities.Cut — TypeCutContains data for a Ticra "Tabulated Pattern Data" object. Note that a single Cut instance contains all of the cuts for a single frequency.
Fields
ncomp::Int: Number of polarization components (2 or 3)icut::Int: 1 for standard constant ϕ polar cuts, 2 for conical, constant θ cuts.icomp::Int: Polarization control parameter. 1 for Eθ and Eφ, 2 for RHCP and LHCP, 3 for Co and Cx (Ludwig 3).text::Vector{String}: Identification text for each constant angle cut.theta::Tt<:AbstractRange: The theta values (in degrees) stored in the cut.phi::Tp<:AbstractRange: The phi values (in degrees) stored in the cut.evec: Matrix of complex field vectors for the two or three polarization components.
TicraUtilities.Exi — TypeExiContains data from a Ticra excitation file.
Fields
header::Vector{String}: Header strings. One element for each line of the header.ampdb::Vector{Float64}: Excitation amplitudes in dB.phsdeg::Vector{Float64}: Excitation phases in degrees.id::Vector{String}: Excitation ID strings.
TicraUtilities.SPHQPartition — TypeSPHQPartitionStruct containing all the data read from a Ticra-compatible Q-type SPH file for one frequency. The qsmns field contains the Q coefficients. It is indexed as qsmns[s,m,n] where s ∈ {1,2}, m ∈ -mmax:mmax, n ∈ {1, ..., mmax}. But not all entries are used. The only possible nonzero entries are where n ∈ {max(1, abs(m)), ..., nmax}. Note that if the coefficients stored in the SPH files are Q', and the cofficients stored in a SPHQPartition instance are Q, then Q = sqrt(8π) * conj(Q'), as discussed in the File Formats section of the Ticra official documentation.
Fields
prgtag::String: Program tag and time stamp, from the program that created the file.idstrg::String: Identification text.nthe::Int: Number of θ-samples over 360°. Must be even and ≥ 4.nphi::Int: Number of ϕ-samples over 360°. Must be ≥ 3.nmax::Int: Maximum value for polar indexnin the SW expansion.1 ≤ nmax ≤ nthe÷2.mmax::Int: Maximum value for azimuthal index|m|in the SW expansion. 0 ≤ mmax ≤ min((nphi-1)÷2, nmax).t4::String: Dummy string.t5::String: Dummy string. If parsed, contains 5 real numbers.t6::String: Dummy string. If parsed, contains 5 real numbers.t7::String: Dummy string.t8::String: Dummy string.qsmns::OffsetArray{ComplexF64, 3}: Contains the Q coefficients, assuming exp(jωt) time variation and using the Ticra normalization convention. The array has axes(1:2, -mmax:mmax, 1:nmax), meaning that some of the entries (those withn < abs(m)) are always zero.powerms::OffsetArray{Float64, 1}: The vector has axes(0:mmax)and themth element contains the total power (one-half the sum of the magnitude squared) of all modes with±mas the m index.
TicraUtilities.Station — TypeStationA struct that contains data for a single partition of a Ticra-compatible optimization station file.
Fields
npoint::Int: The number of stations in the partition.u::Vector{Float64}: A vector of lengthnpointcontaining the $u$ values (unitless direction cosine along $x$) of each of the optimization stations in the partition.v::Vector{Float64}: A vector of lengthnpointcontaining the $v$ values (unitless direction cosine along $y$) of each of the optimization stations in the partition.goal::Vector{Float64}: A vector of lengthnpointcontaining the goal values in dB for each of the optimization stations in the partition.weight::Vector{Float64}: A vector of lengthnpointcontaining the optimization weights for each of the optimization stations in the partition.ipol::Vector{Int}: A vector of lengthnpointcontaining the polarization codes for each of the stations in the partition. The codes have the following meanings:- 1: Linear Ludwig 3 x ("copol") component.
- 2: Linear Ludwig 3 y ("cxpol") component.
- 3: RHCP component.
- 4: LHCP component.
- 5: Major axis of polarization ellipse.
- 6: Minor axis of polarization ellipse.
- 7: The ratio co/cx in dB.
- 8: The ratio cx/co in dB.
- 9: The ratio RHCP/LHCP in dB.
- 10: The ratio LHCP/RHCP in dB.
- 11: The ratio major/minor in dB.
- 12: The ratio minor/major in dB.
- 13: The total power 10*log₁₀(‖E⃗‖²) in dB.
rot::Vector{Float64}: A vector of lengthnpointcontaining the rotation in degrees of linear polarization reference for each of the optimization stations.att::Vector{Float64}: A vector of lengthnpointcontaining the attenuation in dB ≥ 0 relative to the sub-satellite point.ID::Vector{String}: A vector of lengthnpointcontaining the ID string for each of the optimization stations in the partition. May be empty.
TicraUtilities.Surface — TypeSurfaceStruct containing all the data read from a Ticra-compatible regular x-y grid surface file.
Fields
text::String: Descriptive string.idstrg::String: Identification text.x::AbstractRange: The x coordinates at which the surface is sampled.y::AbstractRange: The y coordinates at which the surface is sampled.z::Matrix{Float64}:z[i,j]contains the surface $z$-coordinate sampled at $x$ =x[i]and $y$ =y[j].
TicraUtilities.TEPperiodic — TypeTEPperiodic <: TEPStruct containing the data from a TICRA-compatible Tabulated Electrical Properties (TEP) file for a periodic unit cell. Note that TEP files containing geometrical parameter sweeps are not yet supported.
Fields
name::String: Object name of the periodic unit cell.class::String: Class name of the periodic unit cell.theta: AnAbstractRangecontaining the θ values in units specified by fieldatunit.phi: AnAbstractRangecontaining the ϕ values in units specified by fieldatunit..freqs: A vector of frequencies, each element of which is aUnitfulquantity. The elements may or may not all share the same frequency units.sff::Array{ComplexF64, 5}: Contains reflection coefficients for front surface incidence.sfr::Array{ComplexF64, 5}: Contains transmission coefficients for rear surface incidence.srf::Array{ComplexF64, 5}: Contains transmission coefficients for front surface incidence.srr::Array{ComplexF64, 5}: Contains reflection coefficients for rear surface incidence.
For s representing any of the fields sff, sfr, srf, and srr, size(s) = (2, 2, length(theta), length(phi), length(freqs)), and the 2×2 matrix s[:,:,i,j,k] is arranged in the order [sTETE sTETM; sTMTE sTMTM].
See Also
TEPscatter
TicraUtilities.TEPscatter — TypeTEPscatter <: TEPStruct containing the data for a single frequency from a TICRA-compatible Tabulated Electrical Properties (TEP) file for a scattering surface (i.e., the old-style, original definition of TEP file in use since GRASP8.)
Fields
title::String: Contains the title string.theta: AnAbstractRangecontaining the θ values in degrees.phi: AnAbstractRangecontaining the ϕ values in degrees.sff::Array{ComplexF64, 4}: Contains reflection coefficients for front surface incidence.sfr::Array{ComplexF64, 4}: Contains transmission coefficients for rear surface incidence.srf::Array{ComplexF64, 4}: Contains transmission coefficients for front surface incidence.srr::Array{ComplexF64, 4}: Contains reflection coefficients for rear surface incidence.
For s representing any of the fields sff, sfr, srf, and srr, size(s) = (2, 2, length(theta), length(phi)), and the 2×2 matrix s[:,:,i,j] is arranged in the order [sθθ sθϕ; sϕθ sϕϕ].
See Also
TEPperiodic
TicraUtilities.TorObj — TypeTorObjStruct containing a single, general, TOR object.
Fields
name::String: The name of the TOR object. Example: "my_rim".objtype::String: The type of the TOR object. Example: "tabulatedrimxy".propname::Vector{String}: Names of object properties. Example: ["filename", "unit", "numberofpoints", "translation", "polarorigin"]propval::Vector{String}: Values of object properties corresponding to the names inpropnames. Example: ["h9mscalloped_rim.xyz", "in", "112", "struct(x: 0.0 in, y: 0.0 in)", "struct(status: automatic, x: 0.0 in, y: 0.0 in)"]
Note that the "values" in propval are simply strings that must be parsed for any numeric values that might be present.
See Also
parse_tor_struct
This page was generated using Literate.jl.