PSSFSS Public API Documentation
PSSFSS.analyze — Functionresult = analyze(strata, flist, steering; outlist=[], logfile="pssfss.log", resultfile="pssfss.res", showprogress::Bool=true, fastsweep=true)Analyze a full FSS/PSS structure over a range of frequencies and steering angles/phasings. Generate output files as specified in outlist.
Positional Arguments
strata: A vector ofLayerandRWGSheetobjects. The first and last entries must be of typeLayer.flist: An iterable containing the analysis frequencies in GHz.steering: A length 2NamedTuplecontaining as keys the steering parameter labels and as values the iterables that define the values of steering parameters to be analyzed.one of {
:phi,:ϕ} and one of {:theta,:θ}, orone of {
:psi1,:ψ₁} and one of {:psi2,:ψ₂}.
All steering parameters are input in degrees. Examples of valid
steeringtuples:(θ=0, ϕ=0),(theta=0:10:60, phi=[0, 45]),(theta=40, ϕ=90),(psi1=0, psi2=90),(ψ₁=0, psi2=34.3).The program will analyze while iterating over a triple loop over the two steering parameters and frequency, with frequency in the innermost loop (i.e. varying the fastest). The steering parameter listed first will be in the outermost loop and will therefore vary most slowly.
Keyword Arguments
outlist: A matrix with 2 columns. The first column in each row is a string containing the name of the CSV file to write the output to. The second entry in each row is a tuple generated by the@outputsmacro of theOutputsmodule. The contents of the specified file(s) will be updated as the program completes each analysis frequency.logfile: A string containing the name of the log file to which timing and other information about the run is written. Defaults to"pssfss.log". If this file already exists, it will be overwritten.resultfile: A string containing the name of the results file. Defaults topssfss.res. If this file already exists, it will be overwritten. It is a binary file that contains information (including the generalized scattering matrix) from the analysis performed for each scan condition and frequency. The result file can be post-processed to produce similar or additional outputs that were requested at run time via theoutlistargument.showprogress: If true (default), then show progress bar during execution.fastsweep: If true (default) use an interpolated fast sweep for each frequency loop.
Return Value
result: A vector ofResultobjects, one for each scan angle/frequency combination. This vector can be passed as an input to theextract_resultfunction to obtain any desired performance parameters that are supported by the@outputsmacro.
PSSFSS.Elements.diagstrip — Functiondiagstrip(;P::Real, w::Real, orient::Real, Nl::Int, Nw::Int, units::PSSFSSLength, kwargs...)Return a variable of type RWGSheet that contains the triangulation for a rectangular strip inside a square unit cell, with the strip centerline coincident with one of the diagonals of the cell. The strip runs the full length of the diagonal.
Arguments:
All arguments are keyword arguments which can be entered in any order.
Required arguments:
P: The period (i.e. side length) of the square unit cell. Note that the center-center spacing of the strips isP/√2.w: The width of the strip.orient: The orientation of the strip within the unrotated unit cell in degrees. The only valid values are45for a strip running from lower left to upper right and-45for a strip running from lower right to upper left.units: Length units (mm,cm,inch, ormil)NlandNw: Number of line segments along the length and width of the strip, for dividing up the strip into rectangles, which are triangulated by adding a diagonal to each rectangle. These arguments are actually used for triangulating the central, rectangular portion of the strip. The ends of the strip are tapered in the form of right, isosceles triangles, to conform to the boundaries of the square unit cell. These triangular "end-caps" are triangulated using an unstructured mesh.
Optional arguments:
class::Char='J'Specify the class, either'J'or'M'.. If'J', the unknowns are electric surface currents, as used to model a wire or metallic patch-type FSS. If'M', the unknowns are magnetic surface currents, as used to model a slot or aperture in a perfectly conducting plane.dx::Real=0.0,dy::Real=0.0: These specify the offsets in the x and y directions applied to the entire unit cell and its contents. Length units are as specified in theunitskeyword.rot::Real=0.0: Counterclockwise rotation angle in degrees applied to the entire unit cell and its contents. This rotation is applied prior to any offsets specified indxanddy.Zsheet::Complex=0.0: The frequency-independent surface impedance of the FSS conductor in units of [Ω]. May only be specified for a sheet of class'J'. IfZsheetis specified, thensigma(orσ) may not be specified. )sigmaorσ: DC, bulk conductivity [S/m]. Only allowed for sheets of class'J'. Cannot be simultaneously specified withZsheet. Is used withRqby PSSFSS to calculate an effective sheet surface impedance at each frequency, using the Gradient Model (Grujić 2022).Rq=0.0: RMS surface roughness [m]. Only legal for class'J'. Only used ifsigma(orσ) is also specified. In that case is is used along withsigmato calculate a frequency-dependent sheet impedance using the Gradient Model. The default value of 0 denotes a smooth surface.disttype::Symbol=:normal: Probability distrubution type for surface roughness. defaults to:normal. The other legal value is:rayleigh.fufp::Bool: This keyword is not usually required.fufpis mnemonic for "Find Unique Face Pairs". If true, the code will search the triangulation for classes of triangle pairs that are the equivalent in the toeplitz sense. I.e., if triangle pairs (A,B) and (C,D) belong to the same equivalence class, the six vertices in the pair (A,B) can be made to coincide with those of pair (C,D) by a simple translation. If there are many such equivalent pairs, a significant decrease in matrix fill time ensues by exploiting the equivalence. The tradeoff is the time needed to identify them. The default value istruefor thestrip,diagstrip,meander,manji,loadedcross,jerusalemcross,pixels, 4-sidedpolyringstyles (those employing structured meshes), andsympixels, andfalsefor the remaining styles (those employing unstructured meshes).
PSSFSS.Sheets.edgecount — Functionedgecount(s::RWGSheet)Return the number of triangle edges in the RWGSheet triangulation.
PSSFSS.Sheets.export_sheet — Functionexport_sheet(fname::AbstractString, sheet::RWGSheet, export_type)Export an RWGSheet triangulation to an STL CAD file. export_type may be either STL_ASCII or STL_BINARY.
PSSFSS.Outputs.extract_result — Functionextract_result(r::Result, ops::NTuple{N,Outfun}) --> Row Matrix
extract_result(r::AbstractVector{Result}, ops::NTuple{N,Outfun}) --> MatrixReturn a matrix of outputs extracted from a Result instance or vector. ops is a NTuple as returned by the @outputs macro.
Example
results = analyze(...)
ops = @outputs FGHz s11dB(h,h) s11ang(h,h)
data = extract_result(results, ops)
# or data = extract_result(results[1], ops) # returns a single rowextract_result(fname::AbstractString, ops::Tuple) --> MatrixReturn a matrix of outputs extracted from a results file. ops is a Tuple returned by the @outputs macro.
Example
ops = @outputs FGHz S11DB(H,H) S11ANG(H,H)
data = extract_result("pssfss.res", ops)PSSFSS.Sheets.facecount — Functionfacecount(s::RWGSheet)Return the number of triangle faces in the RWGSheet triangulation.
PSSFSS.Elements.jerusalemcross — Functionjerusalemcross(;P::Real, L1::Real, L2::Real, A::Real, B::Real, w::Real,
ntri::Int, units::PSSFSSLength, kwargs...)Description:
Create a variable of type RWGSheet that contains the triangulation for a "Jerusalem cross" type of geometry. The returned value has fields s₁, s₂, β₁, β₂, ρ, e1, e2, fv, fe, and fr properly initialized.
The following "ascii art" attempts to show the definitions of the geometrical parameters P, L1, L2, A, B, and w. Note that the structure is supposed to be symmetrical wrt reflections about its horizontal and vertical centerlines, and wrt reflections through a line oriented at a 45 degree angle wrt the x-axis.
┏━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┓
┃ ┃ _______
┃ ┌────────────────────────┐ ┃ ↑
┃ │ ┌───────────────────┐ │ ┃ │
┃ │ └───────┐ ┌──────┘ │ ┃ │
┃ └──────┐ │ │ ┌───────┘ ┃ │
┃ │ │ │ │ ┃ │
┃ ┌───────┐ │ │ │ │ ┌──────┐ ┃ │
┃ │ ┌─┐ │ │ │ │ │ │ ┌──┐ │ ┃ │
┃ │ │ │ │ │ │ →│ │← w │ │ │ │ ┃ │
┃ │ │ │ │ │ │ │ │ │ │ │ │ ┃ │
┃ │ │ │ └───────────┘ │ │ └────────────┘ │ │ │ ┃ │
┃ │ │ └─────────────────┘ └────────────────┘ │ │ ┃
┃ │ │ │ │ ┃ L1
┃ │ │ ┌─────────────────┐ ┌────────────────┐ │ │ ┃
┃ │ │ │ ┌───────────┐ │ │ ┌────────────┐ │ │ │ ┃ │
┃ │ │ │ │ │ │ │ │ │ │ │ │ ┃ │
┃ │ │ │ │ │ │ │ │ │ │ │ │ ┃ │
┃ │ └─┘ │ →│ │ │ │← L2 B →│ └──┘ │← ┃ │
┃ └───────┘ │ │ │ │ └──────┘ ┃ │
┃ │ │ │ │ ┃ │
┃ ┌──────┘ │ │ └───────┐ ┃ │
┃ │ ┌───────┘ └──────┐ │ ┃ │
┃ │ └───────────────────┘ │ ┃ │
┃ └────────────────────────┘ ┃ ___↓___
┃ |<───────── A ──────────>| ┃
┗━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━━┛
|<─────────────────────── P ───────────────────────────>|Arguments:
All arguments are keyword arguments which can be entered in any order.
Required arguments:
P: The period, i.e. the side length of the square unit cell.L1,L2,A,B,w: Geometrical parameters as defined above. Note that it is permissible to specifyw ≥ L2/2and/orw ≥ B/2in which case the respective region will be filled in solidly with triangles. If both conditions hold, then the entire structure will be filled in (i.e., singly-connected). In that case theL2andBdimensions will be used for the respective widths of the arms, andwwill not be used.units: Length units (mm,cm,inch, ormil)ntri: The desired total number of triangles. This is a guide/request, the actual number will likely be different.
Optional arguments:
class::Char='J'Specify the class, either'J'or'M'.. If'J', the unknowns are electric surface currents, as used to model a wire or metallic patch-type FSS. If'M', the unknowns are magnetic surface currents, as used to model a slot or aperture in a perfectly conducting plane.dx::Real=0.0,dy::Real=0.0: These specify the offsets in the x and y directions applied to the entire unit cell and its contents. Length units are as specified in theunitskeyword.rot::Real=0.0: Counterclockwise rotation angle in degrees applied to the entire unit cell and its contents. This rotation is applied prior to any offsets specified indxanddy.Zsheet::Complex=0.0: The frequency-independent surface impedance of the FSS conductor in units of [Ω]. May only be specified for a sheet of class'J'. IfZsheetis specified, thensigma(orσ) may not be specified. )sigmaorσ: DC, bulk conductivity [S/m]. Only allowed for sheets of class'J'. Cannot be simultaneously specified withZsheet. Is used withRqby PSSFSS to calculate an effective sheet surface impedance at each frequency, using the Gradient Model (Grujić 2022).Rq=0.0: RMS surface roughness [m]. Only legal for class'J'. Only used ifsigma(orσ) is also specified. In that case is is used along withsigmato calculate a frequency-dependent sheet impedance using the Gradient Model. The default value of 0 denotes a smooth surface.disttype::Symbol=:normal: Probability distrubution type for surface roughness. defaults to:normal. The other legal value is:rayleigh.fufp::Bool: This keyword is not usually required.fufpis mnemonic for "Find Unique Face Pairs". If true, the code will search the triangulation for classes of triangle pairs that are the equivalent in the toeplitz sense. I.e., if triangle pairs (A,B) and (C,D) belong to the same equivalence class, the six vertices in the pair (A,B) can be made to coincide with those of pair (C,D) by a simple translation. If there are many such equivalent pairs, a significant decrease in matrix fill time ensues by exploiting the equivalence. The tradeoff is the time needed to identify them. The default value istruefor thestrip,diagstrip,meander,manji,loadedcross,jerusalemcross,pixels, 4-sidedpolyringstyles (those employing structured meshes), andsympixels, andfalsefor the remaining styles (those employing unstructured meshes).structuredtri::Bool=true: If true, use a structured mesh for the triangulation. If false, the unstructured mesh generator that was standard up to PSSFSS version 1.2 will be used. A structured mesh can be analyzed more efficiently, but the number of triangles created by the unstructured mesh generator is usually closer tontrithan the number for the structured mesh generator.
PSSFSS.Layers.Layer — TypeLayer <: AnyAn instance of the Layer type represents a single dielectric layer of the physical structure. It contains the electrical properties of the dielectric layer. For layers not included in a Gblock, an instance of Layer also specifies the periodicity (via the reciprocal lattice vectors) and stores the mode constants for the Floquet modes present in the layer.
Layer(;width::0mm, ϵᵣ=1.0, tanδ=0.0, μᵣ=1.0, mtanδ=0.0)
Layer(;width::0mm, epsr=1.0, tandel=0.0, mur=1.0, mtandel=0.0)Create a Layer instance with the specified electrical properties. All arguments are optional keyword arguments with default values as shown above. They can be supplied in any order. Typically the first and last Layer in a composite FSS or PSSFSS structure are generated by a call to Layer() to represent the semi-infinite vacuum regions surrounding the structure. Using the default width of 0mm means that the phase reference planes are located at the surfaces just adjacent to these semi-infinite end Layers.
Arguments
width: The layer width (i.e. thickness) expressed as aUnitfullength quantity. For convenience the following unit suffixes are exported by this module:m,cm,mil,inch, so one can specify, e.g.,width=20mil. Note thatwidthcan be negative for the first and/or final layer of the composite structure, which has the effect of shifting the phase reference plane towards the interior of the composite structure. This is sometimes needed when interfacing with other programs, such as TEP file generation for Ticra'sGRASPprogram.ϵᵣorepsr: Relative permittivity of the dielectric.tanδortandel: Loss tangent (electrical) of the dielectric.μᵣormur: Relative permeability of the dielectric.mtanδormtandel: Loss tangent (magnetic) of the dielectric.
PSSFSS.Elements.loadedcross — Functionloadedcross(;s1::Vector{<:Real}, s2::Vector{<:Real}, L1::Real, L2::Real, w::Real,
ntri::Int, units::PSSFSSLength, kwargs...)Description:
Create a variable of type RWGSheet that contains the triangulation for a "loaded cross" type of geometry. The returned value has fields s₁, s₂, β₁, β₂, ρ, e1, e2, fv, fe, and fr properly initialized.
The following (very poor) "ascii art" attempts to show the definitions of the geometrical parameters L1, L2 and w. Note that the structure is supposed to be symmetrical wrt reflections about its horizontal and vertical centerlines, and wrt reflections through a line oriented at a 45 degree angle wrt the x-axis.
^ ----------------
| | _________ |
| | | | |
| | | | |
| | | -->| |<--- W
| | | | |
| | | | |
| ------------ | | -------------
| | |-----------| |------------| |
| | | | |
L1 | | | |
| | | | |
| | | | |
| | ------------ ------------ |
| |----------- | | ------------|
| | | | |
| | | | |
| | | | |
| | | | |
| | |________| |
| | |
V ----------------
<---- L2 ------>Arguments:
All arguments are keyword arguments which can be entered in any order.
Required arguments:
s1ands2: 2-vectors containing the unit cell lattice vectors.L1,L2,w: Geometrical parameters as defined above. Note that it is permissible to specifyw ≥ L2/2in which case a solid (i.e., singly-connected) cross will be generated. In that case theL2dimension will be used for the width of the cross pieces.units: Length units (mm,cm,inch, ormil)ntri: The desired total number of triangles. This is a guide/request, the actual number will likely be different.
Optional arguments:
orient::Real=0.0: Counterclockwise rotation angle in degrees used to locate the initial vertex of the loaded cross. The default is to locate the vertex on the positive x-axis.class::Char='J'Specify the class, either'J'or'M'.. If'J', the unknowns are electric surface currents, as used to model a wire or metallic patch-type FSS. If'M', the unknowns are magnetic surface currents, as used to model a slot or aperture in a perfectly conducting plane.dx::Real=0.0,dy::Real=0.0: These specify the offsets in the x and y directions applied to the entire unit cell and its contents. Length units are as specified in theunitskeyword.rot::Real=0.0: Counterclockwise rotation angle in degrees applied to the entire unit cell and its contents. This rotation is applied prior to any offsets specified indxanddy.Zsheet::Complex=0.0: The frequency-independent surface impedance of the FSS conductor in units of [Ω]. May only be specified for a sheet of class'J'. IfZsheetis specified, thensigma(orσ) may not be specified. )sigmaorσ: DC, bulk conductivity [S/m]. Only allowed for sheets of class'J'. Cannot be simultaneously specified withZsheet. Is used withRqby PSSFSS to calculate an effective sheet surface impedance at each frequency, using the Gradient Model (Grujić 2022).Rq=0.0: RMS surface roughness [m]. Only legal for class'J'. Only used ifsigma(orσ) is also specified. In that case is is used along withsigmato calculate a frequency-dependent sheet impedance using the Gradient Model. The default value of 0 denotes a smooth surface.disttype::Symbol=:normal: Probability distrubution type for surface roughness. defaults to:normal. The other legal value is:rayleigh.fufp::Bool: This keyword is not usually required.fufpis mnemonic for "Find Unique Face Pairs". If true, the code will search the triangulation for classes of triangle pairs that are the equivalent in the toeplitz sense. I.e., if triangle pairs (A,B) and (C,D) belong to the same equivalence class, the six vertices in the pair (A,B) can be made to coincide with those of pair (C,D) by a simple translation. If there are many such equivalent pairs, a significant decrease in matrix fill time ensues by exploiting the equivalence. The tradeoff is the time needed to identify them. The default value istruefor thestrip,diagstrip,meander,manji,loadedcross,jerusalemcross,pixels, 4-sidedpolyringstyles (those employing structured meshes), andsympixels, andfalsefor the remaining styles (those employing unstructured meshes).structuredtri::Bool=true: If true, use a structured mesh for the triangulation. If false, the unstructured mesh generator that was standard up to PSSFSS version 1.2 will be used. A structured mesh can be analyzed more efficiently, but the number of triangles created by the unstructured mesh generator is usually closer tontrithan the number for the structured mesh generator.
PSSFSS.Elements.manji — Functionfunction manji(; L1, L2, L3, w, s1, s2, ntri, units, a=0, w2=0, CCW=false, orient=0, kwargs...)Description:
Create a variable of type RWGSheet that contains the triangulation for a "manji" (Japanese for swastica shape) type of geometry. The returned value has fields s₁, s₂, β₁, β₂, ρ, e1, e2, fv, fe, and fr properly initialized.
Arguments:
All arguments are keyword arguments which can be entered in any order.
Required arguments:
L1,L2,L3,w: Geometrical parameters defined in the diagram at
Note that if a≤wthen the center square shown in the figure will not be present. Similarly, ifL3≤2*wthen the bent portions of the arms will consist of a single strip of widthL3(without any gap in the middle).s1ands2: 2-vectors containing the unit cell lattice vectors.units: Length units (mm,cm,inch, ormil)ntri: The desired total number of triangles. This is a guide/request, the actual number will likely be different.
Optional arguments:
a::Real=0: A geometrical parameter defined in the above referenced diagram. Ifa≤wthen the center square shown in that figure will be absent, and the arms will continue uninterrupted to the center of the structure.w2::Real=0: The width of the square ring border shown in the above-referenced diagram. Ifw2≤ 0 the square ring will not be included in the triangulation. Note thatw2 > 0is only allowed for square unit cells.L4: The outer side length of the square ring border.0 < L4 ≤ norm(s1). If not specified, whenw2 > 0, the default value forL4is the unit cell square dimension. It is the user's responsibility to ensure thatL4is large enough to prevent the square ring from interfering with other parts of themanjistructure.CCW::Bool=false: By default, the chiral structure has a clockwise sense. IfCCWistrue, the structure will be counter-clockwise.orient::Real=0.0: Counterclockwise rotation angle in degrees applied to the structure within the unrotated unit cell. Note that the outer square ring present whenw2 > 0will not be rotated by use of a nonzeroorientvalue.class::Char='J'Specify the class, either'J'or'M'.. If'J', the unknowns are electric surface currents, as used to model a wire or metallic patch-type FSS. If'M', the unknowns are magnetic surface currents, as used to model a slot or aperture in a perfectly conducting plane.dx::Real=0.0,dy::Real=0.0: These specify the offsets in the x and y directions applied to the entire unit cell and its contents. Length units are as specified in theunitskeyword.rot::Real=0.0: Counterclockwise rotation angle in degrees applied to the entire unit cell and its contents. This rotation is applied prior to any offsets specified indxanddy.Zsheet::Complex=0.0: The frequency-independent surface impedance of the FSS conductor in units of [Ω]. May only be specified for a sheet of class'J'. IfZsheetis specified, thensigma(orσ) may not be specified. )sigmaorσ: DC, bulk conductivity [S/m]. Only allowed for sheets of class'J'. Cannot be simultaneously specified withZsheet. Is used withRqby PSSFSS to calculate an effective sheet surface impedance at each frequency, using the Gradient Model (Grujić 2022).Rq=0.0: RMS surface roughness [m]. Only legal for class'J'. Only used ifsigma(orσ) is also specified. In that case is is used along withsigmato calculate a frequency-dependent sheet impedance using the Gradient Model. The default value of 0 denotes a smooth surface.disttype::Symbol=:normal: Probability distrubution type for surface roughness. defaults to:normal. The other legal value is:rayleigh.fufp::Bool: This keyword is not usually required.fufpis mnemonic for "Find Unique Face Pairs". If true, the code will search the triangulation for classes of triangle pairs that are the equivalent in the toeplitz sense. I.e., if triangle pairs (A,B) and (C,D) belong to the same equivalence class, the six vertices in the pair (A,B) can be made to coincide with those of pair (C,D) by a simple translation. If there are many such equivalent pairs, a significant decrease in matrix fill time ensues by exploiting the equivalence. The tradeoff is the time needed to identify them. The default value istruefor thestrip,diagstrip,meander,manji,loadedcross,jerusalemcross,pixels, 4-sidedpolyringstyles (those employing structured meshes), andsympixels, andfalsefor the remaining styles (those employing unstructured meshes).
PSSFSS.Elements.meander — Functionmeander(;a::Real, b::Real, h::Real, w1::Real, w2::Real, ntri::Int,
units::PSSFSSLength, orient=0, kwarg...) --> sheet::RWGSheetDescription:
Return a variable of type RWGSheet that contains the triangulation for a meanderline strip. The returned sheet has the components s₁, s₂, β₁, β₂, ρ, e1, e2, fv, fe, and fr properly initialized. Geometrical parameters are shown in the following diagram:
- - - - - - - - - - - - - - - - - - - - - - - - - ^
| | |
| | |
| | |
| | |
| | |
| <-------- a/2 -------> | |
| (center-to-center) | |
| | |
| ---------------------------- | ^ ^ b
| | | | w2 | |
| | | | | | |
| | ----------------------- | | v | |
| | | | | | |
| -->| |<--w1 w1-->| |<-- | h |
------------ | | ------------ ^ |
| | | | w2 | |
| | | | | | |
------------ - - - - - - - - - - - --------------- v v v
<-------------------- a ------------------------->a and b are unit cell dimensions. w1 and w2 are the widths of the vertical and horizontal strips, resp. h is the total height of the meander.
A nicer diagram: 
Arguments:
All arguments are keyword arguments which can be entered in any order.
Required arguments:
a,b,h,w1,w2: Geometrical parameters as defined above.units: Length units (mm,cm,inch, ormil)ntri: The desired total number of triangles. This is a guide, the actual number will likely be different.
Optional arguments:
orient::Real=0.0: Counterclockwise rotation angle in degrees used to rotate the meanderline orientation within the unrotated unit cell. Nonzero values are allowed only when the unit cell is a square (i.e.a==b). The only allowable values are positive or negative multiples of 90.class::Char='J'Specify the class, either'J'or'M'.. If'J', the unknowns are electric surface currents, as used to model a wire or metallic patch-type FSS. If'M', the unknowns are magnetic surface currents, as used to model a slot or aperture in a perfectly conducting plane.dx::Real=0.0,dy::Real=0.0: These specify the offsets in the x and y directions applied to the entire unit cell and its contents. Length units are as specified in theunitskeyword.rot::Real=0.0: Counterclockwise rotation angle in degrees applied to the entire unit cell and its contents. This rotation is applied prior to any offsets specified indxanddy.Zsheet::Complex=0.0: The frequency-independent surface impedance of the FSS conductor in units of [Ω]. May only be specified for a sheet of class'J'. IfZsheetis specified, thensigma(orσ) may not be specified. )sigmaorσ: DC, bulk conductivity [S/m]. Only allowed for sheets of class'J'. Cannot be simultaneously specified withZsheet. Is used withRqby PSSFSS to calculate an effective sheet surface impedance at each frequency, using the Gradient Model (Grujić 2022).Rq=0.0: RMS surface roughness [m]. Only legal for class'J'. Only used ifsigma(orσ) is also specified. In that case is is used along withsigmato calculate a frequency-dependent sheet impedance using the Gradient Model. The default value of 0 denotes a smooth surface.disttype::Symbol=:normal: Probability distrubution type for surface roughness. defaults to:normal. The other legal value is:rayleigh.fufp::Bool: This keyword is not usually required.fufpis mnemonic for "Find Unique Face Pairs". If true, the code will search the triangulation for classes of triangle pairs that are the equivalent in the toeplitz sense. I.e., if triangle pairs (A,B) and (C,D) belong to the same equivalence class, the six vertices in the pair (A,B) can be made to coincide with those of pair (C,D) by a simple translation. If there are many such equivalent pairs, a significant decrease in matrix fill time ensues by exploiting the equivalence. The tradeoff is the time needed to identify them. The default value istruefor thestrip,diagstrip,meander,manji,loadedcross,jerusalemcross,pixels, 4-sidedpolyringstyles (those employing structured meshes), andsympixels, andfalsefor the remaining styles (those employing unstructured meshes).
PSSFSS.Sheets.nodecount — Functionnodecount(s::RWGSheet)Return the number of unique triangle vertices in the RWGSheet triangulation.
PSSFSS.Outputs.@outputs — Macro@outputs(args...)Convert list of user output requests to a vector of functors that generate the requested outputs when applied to a Result instance. In the conversion process, replace lower case letters with upper case.
Examples
julia> output = @outputs FGHz θ ϕ s11db(te,te) S11ang(Te,te)
julia> output = @outputs FGHz theta phi s21db(R,H) ARdB21(H) ARdB11(v)PSSFSS.Elements.pecsheet — Functionpecsheet()Return a variable of type RWGSheet that contains a perfect electric conducting sheet (i.e. an "E-wall").
PSSFSS.Elements.pixels — Functionpixels(; P, patternmat, units, kwargs...)Description:
Create a variable of type RWGSheet that contains the triangulation for an arbitrarily pixelated square unit cell. The returned value has fields s₁, s₂, β₁, β₂, ρ, e1, e2, fv, fe, and fr properly initialized.
Arguments:
All arguments are keyword arguments which can be entered in any order.
Required arguments:
P::Real > 0: The side length of the square unit cell, specified in units defined by theunitskeyword argument.patternmat::AbstractMatrix{<:Integer}: A square matrix, consisting solely of 1's and 0's. The matrix entries control the metallization pattern in the unit cell, with a 1 ortruevalue denoting a metallized pixel, and a 0 orfalsevalue indicating no metallization. The(i,j)entry corresponds to the pixel centered at(x,y) = ((j-1/2)d, P-(i-1/2)d), whered = P / size(patternmat, 1).units: Length units forP(eithermm,cm,inch, ormil).
Optional arguments:
pdiv::Int = 1: The number of "chops" or subdivisions applied to each square pixel side when forming the geometry triangulation. A value of1(the default) means that the pixels are not subdivided any further, except for a single diagonal across each pixel to form triangles. A value of2would subdivide each pixel into 4 squares. A diagonal edge would be added to each of the resulting squares to form triangles.quad::Bool=false: Iftrue, each subpixel (or pixel, ifpdivis 1) is divided into four triangles by adding two diagonals. Iffalse(the default), only a single diagonal is added to produce two triangles.sym::Bool = false: Atruevalue states that the pixel matrix has vertical and horizontal mirror symmetry, and consists of an even number of rows (and columns). Iftrue, the diagonals added to the squares to form triangles will exhibit the same left-right and up-down mirror symmetry as the collection oftrue(false) pixel locations. Ifsym=falseandquad=false, then all added diagonals across the unit cell will have the same orientation.symhas no effect (i.e. is redundant) ifquad=truesince in that case the two added diagonals already ensure mirror symmetry of the triangulation.class::Char='M'Specify the class, either'J'or'M'. If'J', the unknowns are electric surface currents in the areas corresponding to1values of the pixels. If'M', the unknowns are magnetic surface currents in the areas corresponding to0values of the pixels. It is known that using'J'can result in grossly incorrect results for some geometries where adjacent metallic pixels intersect at only a single point. Therefore, use of only'M'is strongly recommended for mostpixelselements and allsympixelselements.dx::Real=0.0,dy::Real=0.0: These specify the offsets in the x and y directions applied to the entire unit cell and its contents. Length units are as specified in theunitskeyword.rot::Real=0.0: Counterclockwise rotation angle in degrees applied to the entire unit cell and its contents. This rotation is applied prior to any offsets specified indxanddy.Zsheet::Complex=0.0: The frequency-independent surface impedance of the FSS conductor in units of [Ω]. May only be specified for a sheet of class'J'. IfZsheetis specified, thensigma(orσ) may not be specified. )sigmaorσ: DC, bulk conductivity [S/m]. Only allowed for sheets of class'J'. Cannot be simultaneously specified withZsheet. Is used withRqby PSSFSS to calculate an effective sheet surface impedance at each frequency, using the Gradient Model (Grujić 2022).Rq=0.0: RMS surface roughness [m]. Only legal for class'J'. Only used ifsigma(orσ) is also specified. In that case is is used along withsigmato calculate a frequency-dependent sheet impedance using the Gradient Model. The default value of 0 denotes a smooth surface.disttype::Symbol=:normal: Probability distrubution type for surface roughness. defaults to:normal. The other legal value is:rayleigh.fufp::Bool: This keyword is not usually required.fufpis mnemonic for "Find Unique Face Pairs". If true, the code will search the triangulation for classes of triangle pairs that are the equivalent in the toeplitz sense. I.e., if triangle pairs (A,B) and (C,D) belong to the same equivalence class, the six vertices in the pair (A,B) can be made to coincide with those of pair (C,D) by a simple translation. If there are many such equivalent pairs, a significant decrease in matrix fill time ensues by exploiting the equivalence. The tradeoff is the time needed to identify them. The default value istruefor thestrip,diagstrip,meander,manji,loadedcross,jerusalemcross,pixels, 4-sidedpolyringstyles (those employing structured meshes), andsympixels, andfalsefor the remaining styles (those employing unstructured meshes).
PSSFSS.Elements.pmcsheet — Functionpmcsheet()Return a variable of type RWGSheet that contains a perfect magnetic conducting sheet (i.e. an "H-wall").
PSSFSS.Elements.polyring — Functionpolyring(;s1::Vector, s2::Vector, a::Vector, b::Vector, sides::Int ,ntri::Int ,orient::Real, units::PSSFSSLength, kwargs...) --> RWGSheetReturn a variable of type RWGSheet that contains the triangulation for one or more concentric annular regions bounded by polygons.
Arguments:
All arguments are keyword arguments which can be entered in any order.
Required arguments:
units: Length units (mm,cm,inch, ormil)s1ands2: 2-vectors containing the unit cell lattice vectors.aandb: n-vectors (n>=1) of the same length providing the inner and outer radii, respectively of the polygonal rings. Entries inaandbmust be strictly increasing, except for possiblyb[end]as discussed below.b[i] > a[i]∀i ∈ 1:n, except possiblyb[end]as discussed below.a[1]may be zero to denote a solid (non-annular) polygon as the first "ring". It is possible to let the outermost ring to extend completely to the unit cell boundary. This is specified by settingb[end]< 0, in which case for unstructured meshes,-b[end]is interpreted as the number of edges along the shorter of thes1ands2lattice vectors.sides: The number (>= 3) of polygon sides.ntri: The desired total number of triangles distributed among all the annular regions. This is a guide, the actual number will likely be different.
Optional arguments:
orient::Real=0.0: Counterclockwise rotation angle in degrees used to locate the initial vertex of the polygonal rings. The default is to locate the vertex on the positive x-axis.structuredtri::Bool: Defaults totruewhensides==4and false otherwise. Atruevalue is only allowed whensides==4ands1⟂s2. If true, use a structured mesh for the triangulation. If false, the unstructured mesh generator that was standard up to PSSFSS version 1.2 will be used. A structured mesh can be analyzed more efficiently, but the number of triangles created by the unstructured mesh generator is usually closer tontrithan the number for the structured mesh generator.class::Char='J'Specify the class, either'J'or'M'.. If'J', the unknowns are electric surface currents, as used to model a wire or metallic patch-type FSS. If'M', the unknowns are magnetic surface currents, as used to model a slot or aperture in a perfectly conducting plane.dx::Real=0.0,dy::Real=0.0: These specify the offsets in the x and y directions applied to the entire unit cell and its contents. Length units are as specified in theunitskeyword.rot::Real=0.0: Counterclockwise rotation angle in degrees applied to the entire unit cell and its contents. This rotation is applied prior to any offsets specified indxanddy.Zsheet::Complex=0.0: The frequency-independent surface impedance of the FSS conductor in units of [Ω]. May only be specified for a sheet of class'J'. IfZsheetis specified, thensigma(orσ) may not be specified. )sigmaorσ: DC, bulk conductivity [S/m]. Only allowed for sheets of class'J'. Cannot be simultaneously specified withZsheet. Is used withRqby PSSFSS to calculate an effective sheet surface impedance at each frequency, using the Gradient Model (Grujić 2022).Rq=0.0: RMS surface roughness [m]. Only legal for class'J'. Only used ifsigma(orσ) is also specified. In that case is is used along withsigmato calculate a frequency-dependent sheet impedance using the Gradient Model. The default value of 0 denotes a smooth surface.disttype::Symbol=:normal: Probability distrubution type for surface roughness. defaults to:normal. The other legal value is:rayleigh.fufp::Bool: This keyword is not usually required.fufpis mnemonic for "Find Unique Face Pairs". If true, the code will search the triangulation for classes of triangle pairs that are the equivalent in the toeplitz sense. I.e., if triangle pairs (A,B) and (C,D) belong to the same equivalence class, the six vertices in the pair (A,B) can be made to coincide with those of pair (C,D) by a simple translation. If there are many such equivalent pairs, a significant decrease in matrix fill time ensues by exploiting the equivalence. The tradeoff is the time needed to identify them. The default value istruefor thestrip,diagstrip,meander,manji,loadedcross,jerusalemcross,pixels, 4-sidedpolyringstyles (those employing structured meshes), andsympixels, andfalsefor the remaining styles (those employing unstructured meshes).
PSSFSS.Sheets.read_sheet_data — Functionread_sheet_data(filename::AbstractString)::RWGSheetRead the sheet triangulation and unit cell data from a JLD2 file named in filename.
PSSFSS.Elements.rectstrip — Functionrectstrip(;Lx::Real, Ly::Real, Nx::Int, Ny::Int, Px::Real, Py::Real, units::PSSFSSLength, kwargs...)Return a variable of type RWGSheet that contains the triangulation for a rectangular strip.
Arguments:
All arguments are keyword arguments which can be entered in any order.
Required arguments:
units: Length units (mm,cm,inch, ormil)LxandLy: Lengths of the strip in the x and y directions.PxandPy: Lengths (periods) of the rectangular unit cell in the x and y directions.NxandNy: Number of line segments in the x and y directions, for dividing up the strip into rectangles, which are triangulated by adding a diagonal to each rectangle.
Optional arguments:
orient::Real=0.0: Counterclockwise rotation angle in degrees applied to the strip within the unrotated unit cell. This rotation is applied prior to any offsets specified indxanddyand any unit cell rotation specified byrot.class::Char='J'Specify the class, either'J'or'M'.. If'J', the unknowns are electric surface currents, as used to model a wire or metallic patch-type FSS. If'M', the unknowns are magnetic surface currents, as used to model a slot or aperture in a perfectly conducting plane.dx::Real=0.0,dy::Real=0.0: These specify the offsets in the x and y directions applied to the entire unit cell and its contents. Length units are as specified in theunitskeyword.rot::Real=0.0: Counterclockwise rotation angle in degrees applied to the entire unit cell and its contents. This rotation is applied prior to any offsets specified indxanddy.Zsheet::Complex=0.0: The frequency-independent surface impedance of the FSS conductor in units of [Ω]. May only be specified for a sheet of class'J'. IfZsheetis specified, thensigma(orσ) may not be specified. )sigmaorσ: DC, bulk conductivity [S/m]. Only allowed for sheets of class'J'. Cannot be simultaneously specified withZsheet. Is used withRqby PSSFSS to calculate an effective sheet surface impedance at each frequency, using the Gradient Model (Grujić 2022).Rq=0.0: RMS surface roughness [m]. Only legal for class'J'. Only used ifsigma(orσ) is also specified. In that case is is used along withsigmato calculate a frequency-dependent sheet impedance using the Gradient Model. The default value of 0 denotes a smooth surface.disttype::Symbol=:normal: Probability distrubution type for surface roughness. defaults to:normal. The other legal value is:rayleigh.fufp::Bool: This keyword is not usually required.fufpis mnemonic for "Find Unique Face Pairs". If true, the code will search the triangulation for classes of triangle pairs that are the equivalent in the toeplitz sense. I.e., if triangle pairs (A,B) and (C,D) belong to the same equivalence class, the six vertices in the pair (A,B) can be made to coincide with those of pair (C,D) by a simple translation. If there are many such equivalent pairs, a significant decrease in matrix fill time ensues by exploiting the equivalence. The tradeoff is the time needed to identify them. The default value istruefor thestrip,diagstrip,meander,manji,loadedcross,jerusalemcross,pixels, 4-sidedpolyringstyles (those employing structured meshes), andsympixels, andfalsefor the remaining styles (those employing unstructured meshes).
PSSFSS.Outputs.res2fresnel — Functionres2fresnel(results::Vector{Result}, fresnelfile::AbstractString)
res2fresnel(resultfile::AbstractString, fresnelfile::AbstractString)Create an HFSS-compatible "Fresnel table" file from results, the vector of Result objects returned by the analyze function. If the first positional argument is an AbstractString, it is assumed to be the name of a PSSFSS results file, from which the vector of results will be read.
Since Fresnel tables contain data for only a single ϕ value, if the input result vector contains data for multiple ϕ values, only the value with minimum magnitude will be used.
Fresnel tables may be formatted to contain only reflection coefficients (for a so-called "opaque" structure), or they may contain both reflection and transmission coefficients (a "non-opaque" structure). An opaque structure is one for which the s21 partition of the generalized scattering matrix is identically zero for all frequencies and scan angles. The correct format to be written will be selected automatically by res2fresnel.
Requirements for Fresnel Table Compatibility
The data in results must satisfy the following requirements:
- Incidence angles rather than incremental phasings must be used.
- θ angles must begin at 0 and be uniformly spaced up to the maximum θ value present.
- The increment in θ values must divide evenly into 90.
- If multiple frequencies are present, then they must have a uniform spacing.
A Fresnel table must contain θ values equally spaced between 0 and 90, inclusive. If the results vector provided as input does not contain θ values all the way to 90, then the scattering matrix values corresponding to the maximum provided θ value will be copied into the remaining angular "slots" as necessary to provide a complete Fresnel table.
There are some limitations on the type of unit cell geometry that should be used for creating Fresnel tables. First, a Fresnel table contains data for only a single ϕ value. This means that the geometry being analyzed must be such that the scattering matrix of the structure is essentially independent of ϕ. As a counterexample, a strip grid is not a suitable structure, since its scattering properties are strongly dependent on ϕ. Second, a Fresnel table records only co-polarized (TE → TE and TM → TM) transmission and reflection coefficients. This means that the structure being analyzed must not generate cross-polarized (TE → TM or TM → TE) transmission or reflection coefficients of significant amplitude.
Fresnel tables consider only incidence from a single "front" region. When creating the Fresnel table, the front region is taken to be Region 1 of the PSSFSS model (i.e. the first layer present in the PSSFSS strata vector).
Additional Requirements for Non-Opaque Structures
When used in an HFSS SBR+ model, the scattering properties read from the Fresnel table are applied to a zero-thickness surface, so that the transmitted ray is launched from the same "hit" point of the surface that was encountered by the incident ray. Because of this, the phase reference plane for both input and output ports of the PSSFSS model should be located at this front surface (i.e. the first interface plane in the strata vector). This is accomplished by specifying zero width for the first Layer object (i.e. using Layer() for the first layer), and then specifying the final layer's width to be the negative of the sum of all the other layer widths in the strata vector. The negative width value shifts the output port reference plane to coincide with that of the input port. As an example:
strata = [Layer(), Layer(width=2mm, ϵᵣ=2.2) Layer(width=3.3mm, ϵᵣ=3.0), Layer(width=2mm, ϵᵣ=2.2), Layer(width=-7.3mm)]PSSFSS.Outputs.res2tep — Functionres2tep(results::Vector{Result}; name="tep", class="res2tep") -> t::TEPperiodic
res2tep(resultfile::AbstractString; name="tep", class="res2tep") -> t::TEPperiodic
res2tep(results::Vector{Result}, tepfile::AbstractString; name="tep", class="res2tep") -> t::TEPperiodic
res2tep(resultfile::AbstractString, tepfile::AbstractString; name="tep", class="res2tep") -> t::TEPperiodicConvert a vector of Result elements into a TEPperiodic object, as defined in the TicraUtilities package. If positional argument tepfile is provided, the TEPperiodic object will be saved to this file name as a TICRA-compatible TEP (tabulated electrical properties) file. If the first positional argument is an AbstractString, it is assumed to be the name of a PSSFSS results file, from which the vector of results will be read. The keyword arguments are used to provide values for the same-named fields in the TEP structure.
Requirements for TEP File Compatibility
Because a TEP file contains all of the information of the full 4×4 scattering matrix computed by PSSFSS, there are no limitations on the type of unit cell geometry that can be used for creating TEP files.
TEP files use the concept of "front" and "rear" incidence. When converting a PSSFSS analysis result to TEP format, Region 1 (the first layer in the strata vector) is taken as the "front" incidence region, and Region n (the last layer) is taken to be the "rear" region. Both of these layers should have zero width and assume vacuum electrical parameters. I.e., they should be specified as Layer() in the strata stackup.
results (or resultfile) must contain the results of a PSSFSS analysis sweep over θ and ϕ (and possibly frequency) such that
- Incidence angles θ and ϕ rather than incremental phasings ψ₁ and ψ₂ were used.
- If more than one ϕ value is used, then all ϕ values in the range
0:Δϕ:(360-Δϕ)must be present. - The entire 3-dimensional Cartesian product of all angles and frequencies must be present.
PSSFSS.Sheets.RWGSheet — TypeRWGSheetA type that represents a zero-thickness sheet of periodically patterned metalization. Particular instances are created by calling a constructor function for a specific type of sheet geometry. These include: diagstrip, jerusalemcross, loadedcross, manji, meander, pecsheet, pmcsheet, polyring, rectstrip, sinuous, and splitring.
PSSFSS.Elements.sinuous — Functionsinuous(; arms, b, w, g, sides, ntri, units, s1, s2, kwargs...) --> RWGSheetReturn a variable of type RWGSheet representing a sinuous element as shown in this diagram: 
Arguments:
All arguments are keyword arguments which can be entered in any order.
Required arguments:
arms::Int: The number of arms in the structure.rc::Real > 0: The radius of the central circle.rcmust be greater than or equal tow.b: n-vector (n ≥ 1) providing the outer radii of the polygonal rings. Entries must be positive and strictly increasing, with the difference between adjacent rings exceedingw.w: The width of the traces in the arms.g: A scalar containing the rectangular gap width separating adjacent arms.sides::Int: The number (>= 4) of polygon sides for the background regular annular polygon(s) from which the ring sections are created.ntri::Int: The desired total number of triangles. This is a guide/request, the actual number will likely be different.units: Length units (mm,cm,inch, ormil)s1ands2: 2-vectors containing the unit cell lattice vectors.
Optional arguments:
orient::Real=0.0: Counterclockwise rotation angle in degrees used to locate center of the first arm.w2::Real=0.0: The trace width of the enclosing square loop "rim". Note thatw2 > 0is only permitted for a square unit cell.L2: The outer dimension (i.e. the full side length) of the square "rim" present whenw2 > 0. The user is responsible for choosingL2large enough that the rim does not intefere with the sinuous arms of the structure.L2must be less than or equal to the square unit cell dimension. It defaults to the unit cell dimension if it is not specified.c2::Real=0.0: The outer dimension of the small squares shown in the corners of the enclosing square loop "rim". Ifc2==0then the squares are not included, and the outer loop is a simple square loop.class::Char='J'Specify the class, either'J'or'M'.. If'J', the unknowns are electric surface currents, as used to model a wire or metallic patch-type FSS. If'M', the unknowns are magnetic surface currents, as used to model a slot or aperture in a perfectly conducting plane.dx::Real=0.0,dy::Real=0.0: These specify the offsets in the x and y directions applied to the entire unit cell and its contents. Length units are as specified in theunitskeyword.rot::Real=0.0: Counterclockwise rotation angle in degrees applied to the entire unit cell and its contents. This rotation is applied prior to any offsets specified indxanddy.Zsheet::Complex=0.0: The frequency-independent surface impedance of the FSS conductor in units of [Ω]. May only be specified for a sheet of class'J'. IfZsheetis specified, thensigma(orσ) may not be specified. )sigmaorσ: DC, bulk conductivity [S/m]. Only allowed for sheets of class'J'. Cannot be simultaneously specified withZsheet. Is used withRqby PSSFSS to calculate an effective sheet surface impedance at each frequency, using the Gradient Model (Grujić 2022).Rq=0.0: RMS surface roughness [m]. Only legal for class'J'. Only used ifsigma(orσ) is also specified. In that case is is used along withsigmato calculate a frequency-dependent sheet impedance using the Gradient Model. The default value of 0 denotes a smooth surface.disttype::Symbol=:normal: Probability distrubution type for surface roughness. defaults to:normal. The other legal value is:rayleigh.fufp::Bool: This keyword is not usually required.fufpis mnemonic for "Find Unique Face Pairs". If true, the code will search the triangulation for classes of triangle pairs that are the equivalent in the toeplitz sense. I.e., if triangle pairs (A,B) and (C,D) belong to the same equivalence class, the six vertices in the pair (A,B) can be made to coincide with those of pair (C,D) by a simple translation. If there are many such equivalent pairs, a significant decrease in matrix fill time ensues by exploiting the equivalence. The tradeoff is the time needed to identify them. The default value istruefor thestrip,diagstrip,meander,manji,loadedcross,jerusalemcross,pixels, 4-sidedpolyringstyles (those employing structured meshes), andsympixels, andfalsefor the remaining styles (those employing unstructured meshes).
PSSFSS.Elements.splitring — Functionsplitring(; s1, s2, a, b, sides, ntri, gapwidth, gapcenter, gapangle, units, kwargs...) --> RWGSheetReturn a variable of type RWGSheet similar to a polyring but with zero or more gaps in each concentric annular region.
Arguments:
All arguments are keyword arguments which can be entered in any order.
Required arguments:
units: Length units (mm,cm,inch, ormil)s1ands2: 2-vectors containing the unit cell lattice vectors.aandb: n-vectors (n>=1) of the same length providing the inner and outer radii, respectively of the polygonal rings. Entries inaandbmust be positive and strictly increasing.b[i] > a[i]∀i ∈ 1:n.sides: The number (>= 4) of polygon sides for the background regular annular polygon(s) from which the gaps are removed.gapcenter: A scalar or vector of angles in degrees that define the gap center angular location(s), measured counterclockwise. A scalar implies that all rings have a gap in that same angular location. If a vector, then it must have the same length asaandb, withgapcenter[m]denoting the gap center location for themth ring. Howevergapcenter[m]can be either a scalar (denoting a single gap) or an n-tuple (denoting n gaps in themth ring).gapwidth: A scalar or a vector of the same length asaandbcontaining the gap width(s) for each ring. A width of zero implies that the ring is not split (i.e. there is no gap). If thegapwidthof all rings is zero, then the resulting geometry is similar to apolyring. If a ring is to have multiple gaps, then the widths of the gaps for that ring should be passed as a tuple. For example, suppose there are three rings and the second ring has 2 gaps, with the others having a single gap. Thengapwidth = [0.5, (0.4, 0.6), 0.3]would be an appropriately formatted input in this case. Whengapwidthis specified, the gaps are implemented as if a rectangular region is removed from the annular polygonal rings. Note that only one ofgapwidthandgapanglecan be specified.gapangle: A scalar or vector of the same length asaandbcontaining the angular widths of the gaps in degrees. As withgapwidth, for any rings with multiple gaps, the corresponding entry ingapangleshould be a tuple of the same length as the number of gaps for that ring. Whengapangleis specified, the gap(s) in themth ring is/are formed as if pie-shaped wedge(s) with wedge angle(s)gapangle[m], are removed from the ring(s). The locations and sizes of the tuples ingapanglemust agree with those ingapcenter. Note that only one ofgapangleandgapwidthcan be specified.ntri: The desired total number of triangles distributed among all the annular regions. This is a guide, the actual number will likely be different.
Optional arguments:
orient::Real=0.0: Counterclockwise rotation angle in degrees used to locate the initial vertex of the polygonal rings. The default is to locate the vertex on the ray from the center parallel to the positive x-axis.class::Char='J'Specify the class, either'J'or'M'.. If'J', the unknowns are electric surface currents, as used to model a wire or metallic patch-type FSS. If'M', the unknowns are magnetic surface currents, as used to model a slot or aperture in a perfectly conducting plane.dx::Real=0.0,dy::Real=0.0: These specify the offsets in the x and y directions applied to the entire unit cell and its contents. Length units are as specified in theunitskeyword.rot::Real=0.0: Counterclockwise rotation angle in degrees applied to the entire unit cell and its contents. This rotation is applied prior to any offsets specified indxanddy.Zsheet::Complex=0.0: The frequency-independent surface impedance of the FSS conductor in units of [Ω]. May only be specified for a sheet of class'J'. IfZsheetis specified, thensigma(orσ) may not be specified. )sigmaorσ: DC, bulk conductivity [S/m]. Only allowed for sheets of class'J'. Cannot be simultaneously specified withZsheet. Is used withRqby PSSFSS to calculate an effective sheet surface impedance at each frequency, using the Gradient Model (Grujić 2022).Rq=0.0: RMS surface roughness [m]. Only legal for class'J'. Only used ifsigma(orσ) is also specified. In that case is is used along withsigmato calculate a frequency-dependent sheet impedance using the Gradient Model. The default value of 0 denotes a smooth surface.disttype::Symbol=:normal: Probability distrubution type for surface roughness. defaults to:normal. The other legal value is:rayleigh.fufp::Bool: This keyword is not usually required.fufpis mnemonic for "Find Unique Face Pairs". If true, the code will search the triangulation for classes of triangle pairs that are the equivalent in the toeplitz sense. I.e., if triangle pairs (A,B) and (C,D) belong to the same equivalence class, the six vertices in the pair (A,B) can be made to coincide with those of pair (C,D) by a simple translation. If there are many such equivalent pairs, a significant decrease in matrix fill time ensues by exploiting the equivalence. The tradeoff is the time needed to identify them. The default value istruefor thestrip,diagstrip,meander,manji,loadedcross,jerusalemcross,pixels, 4-sidedpolyringstyles (those employing structured meshes), andsympixels, andfalsefor the remaining styles (those employing unstructured meshes).
PSSFSS.Elements.sympixels — Functionfunction sympixels(; P, nrim, halfnint, patternvec, units, kwargs...)Description:
Create a variable of type RWGSheet that contains the triangulation for a symmetrically pixelated square unit cell. The pattern of metallic pixels has C4 symmetry, as well as left-right and up-down mirror symmetry, as well as each quadrant exhibiting antidiagonal mirror symmetry. The returned value has fields s₁, s₂, β₁, β₂, ρ, e1, e2, fv, fe, and fr properly initialized. The pixels included in the triangulation are determined by the patternvec input vector as described below.
Arguments:
All arguments are keyword arguments which can be entered in any order.
Required arguments:
P::Real > 0: The side length of the square unit cell, specified in units defined by theunitskeyword argument.nrim::Integer: The width of the solid metallic rim placed just inside the unit cell boundary, in units of pixels. A value of0implies that there is no rim.halfnint: Half the number of pixels in the side length for the interior (strictly inside the rim) square pixelated region of the unit cell.patternvec::AbstractVector{<:Integer}: A vector of lengthhalfnint*(halfnint+1)÷2, consisting solely of 1's and 0's. The elements of this vector are mapped to pixels in the irreducible zone of the unit cell as shown in the diagram at
. Within the irreducible zone, pixels corresponding to a value of 1(ortrue) are taken to be areas of metallization, while0orfalsevalues are metal-free (void) areas. This holds for eitherJorMas theclassvalue (see theclassargument below for important limitations).units: Length units forP(eithermm,cm,inch, ormil).
Optional arguments:
pdiv::Int = 1: The number of "chops" or subdivisions applied to each square pixel side when forming the triangulation. A value of1(the default) means that the pixels included in the triangulation (1ortruevalues forclass='J',0orfalsevalues forclass='M') are not subdivided any further, except for a single diagonal across each square pixel to form triangles. A value ofn>1means that each square pixel is first divided inton×nsquare subpixels, after which a single diagonal edge (ifquad=false) is added to each subpixel to form triangles.quad::Bool=false: Iftrue, each subpixel (or pixel, ifpdivis 1) is divided into four triangles by adding two diagonals. Iffalse(the default), only a single diagonal is added to each square to produce two triangles.sym::Bool = false: If true, the diagonals added to the squares will exhibit the same left-right and up-down mirror symmetry as the collection oftrue(false) pixel locations. Iffalseandquad=false, then all added diagonals across the unit cell will have the same orientation.symhas no effect (i.e. is redundant) ifquad=truesince in that case the two added diagonals already ensure mirror symmetry of the triangulation.class::Char='M'Specify the class, either'J'or'M'. If'J', the unknowns are electric surface currents in the areas corresponding to1values of the pixels. If'M', the unknowns are magnetic surface currents in the areas corresponding to0values of the pixels. It is known that using'J'can result in grossly incorrect results for some geometries where adjacent metallic pixels intersect at only a single point. Therefore, use of only'M'is strongly recommended for mostpixelselements and allsympixelselements.dx::Real=0.0,dy::Real=0.0: These specify the offsets in the x and y directions applied to the entire unit cell and its contents. Length units are as specified in theunitskeyword.rot::Real=0.0: Counterclockwise rotation angle in degrees applied to the entire unit cell and its contents. This rotation is applied prior to any offsets specified indxanddy.Zsheet::Complex=0.0: The frequency-independent surface impedance of the FSS conductor in units of [Ω]. May only be specified for a sheet of class'J'. IfZsheetis specified, thensigma(orσ) may not be specified. )sigmaorσ: DC, bulk conductivity [S/m]. Only allowed for sheets of class'J'. Cannot be simultaneously specified withZsheet. Is used withRqby PSSFSS to calculate an effective sheet surface impedance at each frequency, using the Gradient Model (Grujić 2022).Rq=0.0: RMS surface roughness [m]. Only legal for class'J'. Only used ifsigma(orσ) is also specified. In that case is is used along withsigmato calculate a frequency-dependent sheet impedance using the Gradient Model. The default value of 0 denotes a smooth surface.disttype::Symbol=:normal: Probability distrubution type for surface roughness. defaults to:normal. The other legal value is:rayleigh.fufp::Bool: This keyword is not usually required.fufpis mnemonic for "Find Unique Face Pairs". If true, the code will search the triangulation for classes of triangle pairs that are the equivalent in the toeplitz sense. I.e., if triangle pairs (A,B) and (C,D) belong to the same equivalence class, the six vertices in the pair (A,B) can be made to coincide with those of pair (C,D) by a simple translation. If there are many such equivalent pairs, a significant decrease in matrix fill time ensues by exploiting the equivalence. The tradeoff is the time needed to identify them. The default value istruefor thestrip,diagstrip,meander,manji,loadedcross,jerusalemcross,pixels, 4-sidedpolyringstyles (those employing structured meshes), andsympixels, andfalsefor the remaining styles (those employing unstructured meshes).
PSSFSS.Sheets.write_sheet_data — Functionwrite_sheet_data(filename::AbstractString, sheet::RWGSheet)Write the sheet triangulation and unit cell data to a JLD2 file named in filename.