Geometry API

Aircraft Geometry API

AeroFuse.AircraftGeometry.FoilType
Foil(x, y, name = "")
Foil(coordinates, name = "")

Structure consisting of foil coordinates in 2 dimensions with an optional name.

The coordinates should be provided in counter-clockwise format, viz. from the trailing edge of the upper surface to the trailing edge of the lower surface.

source
AeroFuse.AircraftGeometry.HyperEllipseFuselageMethod
HyperEllipseFuselage(;

Define a fuselage based on the following hyperelliptical-cylindrical parameterization.

Nose: Hyperellipse $z(ξ) = (1 - ξ^a)^{(1/a)}$

Cabin: Cylindrical $z(ξ) = 1$

Rear: Hyperellipse $z(ξ) = (1 - ξ^b)^{(1/b)}$

Arguments

  • radius :: Real = 1.: Fuselage radius (m)
  • length :: Real = 6.: Fuselage length (m)
  • x_a :: Real = 1.: Location of front of cabin as a ratio of fuselage length ∈ [0,1]
  • x_b :: Real = 0.: Location of rear of cabin as a ratio of fuselage length ∈ [0,1]
  • c_nose :: Real = 0.: Curvature of nose in terms of hyperellipse parameter $a$
  • c_rear :: Real = 1.: Curvature of rear in terms of hyperellipse parameter $b$
  • d_nose :: Real = 1.: "Droop" or "rise" of nose from front of cabin centerline (m)
  • d_rear :: Real = 1.: "Droop" or "rise" of rear from rear of cabin centerline (m)
  • position :: Vector{Real} = zeros(3): Position (m)
  • angle :: Real = 0.: Angle of rotation (degrees)
  • axis :: Vector{Real} = [0, 1 ,0]: Axis of rotation, y-axis by default
  • affine :: AffineMap = AffineMap(AngleAxis(deg2rad(angle), axis...), position): Affine mapping for the position and orientation via CoordinateTransformations.jl (overrides angle and axis if specified)
source
AeroFuse.AircraftGeometry.WingMethod
Wing(
    chords, 
    foils :: Vector{Foil}, 
    twists, 
    spans, 
    dihedrals, 
    sweeps,
    symmetry = false,
    flip     = false,
    position = zeros(3),
    angle    = 0.,
    axis     = [0.,1.,0.],
)

Definition for a Wing consisting of $N+1$ Foils, their associated chord lengths $c$ and twist angles $ι$, for $N$ sections with span lengths $b$, dihedrals $δ$ and leading-edge sweep angles $Λ_{LE}$, with all angles in degrees. Optionally, specify translation and a rotation in angle-axis representation for defining coordinates in a global axis system. Additionally, specify Boolean arguments for symmetry or reflecting in the $x$-$z$ plane.

Arguments

  • chords :: Vector{Real}: Chord lengths (m)
  • foils :: Vector{Foil} = fill(naca4(0,0,1,2), length(chords)): Foil shapes, default is NACA 0012.
  • spans :: Vector{Real} = ones(length(chords) - 1) / (length(chords) - 1): Span lengths (m), default yields total span length 1.
  • dihedrals :: Vector{Real} = zero(spans): Dihedral angles (deg), default is zero.
  • sweeps :: Vector{Real} = zero(spans): Sweep angles (deg), default is zero.
  • w_sweep :: Real = 0.: Chord ratio for sweep angle e.g., 0 = Leading-edge sweep, 1 = Trailing-edge sweep, 0.25 = Quarter-chord sweep
  • symmetry :: Bool = false: Symmetric in the $x$-$z$ plane
  • flip :: Bool = false: Flip coordinates in the $x$-$z$ plane
  • position :: Vector{Real} = zeros(3): Position (m)
  • angle :: Real = 0.: Angle of rotation (degrees)
  • axis :: Vector{Real} = [0.,1.,0.]: Axis of rotation
  • affine :: AffineMap = AffineMap(AngleAxis(deg2rad(angle), axis...), position): Affine mapping for the position and orientation via CoordinateTransformations.jl (overrides angle and axis if specified)
source
AeroFuse.AircraftGeometry.WingMeshType
WingMesh(
    wing :: AbstractWing, 
    n_span :: Vector{Integer}, n_chord :: Integer;
    span_spacing :: AbstractSpacing = symmetric_spacing(wing)
)

Define a container to generate meshes and panels for a given Wing with a specified distribution of number of spanwise panels, and a number of chordwise panels.

Optionally a combination of AbstractSpacing types (Sine(), Cosine(), Uniform()) can be provided to the named argument span_spacing, either as a singleton or as a vector with length equal to the number of spanwise sections. By default, the combination is [Sine(), Cosine(), ..., Cosine()].

For surface coordinates, the wing mesh will have (nchord - 1) * 2 chordwise panels from TE-LE-TE and (nspan * 2) spanwise panels.

source
AeroFuse.AircraftGeometry.WingSectionMethod
WingSection(; 
    area, aspect, taper
    dihedral, sweep, w_sweep,
    root_twist, tip_twist,
    position, angle, axis,
    symmetry, flip,
    root_foil, tip_foil,
    root_control, tip_control,
)

Define a Wing in the $x$-$z$ plane, with optional Boolean arguments for symmetry and flipping in the plane.

Arguments

  • area :: Real = 1.: Area (m²)
  • aspect :: Real = 6.: Aspect ratio
  • taper :: Real = 1.: Taper ratio of tip to root chord
  • dihedral :: Real = 1.: Dihedral angle (degrees)
  • sweep :: Real = 0.: Sweep angle (degrees)
  • w_sweep :: Real = 0.: Chord ratio for sweep angle e.g., 0 = Leading-edge sweep, 1 = Trailing-edge sweep, 0.25 = Quarter-chord sweep
  • root_twist :: Real = 0.: Twist angle at root (degrees)
  • tip_twist :: Real = 0.: Twist angle at tip (degrees)
  • root_foil :: Foil = naca4((0,0,1,2)): Foil at root
  • tip_foil :: Foil = root_foil: Foil at tip. Defaults to root foil.
  • root_control :: NTuple{2} = (0., 0.75): (Angle, hinge ratio) for adding a control surface at the root.
  • tip_control :: NTuple{2} = root_control: (Angle, hinge ratio) for adding a control surface at the tip. Defaults to root control's settings.
  • symmetry :: Bool = false: Symmetric in the $x$-$z$ plane
  • flip :: Bool = false: Flip coordinates in the $x$-$z$ plane
  • position :: Vector{Real} = zeros(3): Position (m)
  • angle :: Real = 0.: Angle of rotation (degrees)
  • axis :: Vector{Real} = [0.,1.,0.]: Axis of rotation
  • affine :: AffineMap = AffineMap(AngleAxis(deg2rad(angle), axis...), position): Affine mapping for the position and orientation via CoordinateTransformations.jl (overrides angle and axis if specified)
source
AeroFuse.AircraftGeometry.affineMethod
affine(
    foil :: Foil, 
    angle, vector
)

Perform an affine transformation on the coordinates of a Foil by a 2-dimensional vector $\mathbf v$ and angle $θ$.

source
AeroFuse.AircraftGeometry.camber_CSTFunction
camber_CST(α_c, α_t,
           Δz :: Real,
           n :: Integer = 40)

Define a cosine-spaced foil with $2n$ points using the Class Shape Transformation method on a Bernstein polynomial basis for the camber and thickness coordinates.

The foil is defined by arrays of coefficients $(α_c,~ α_t)$ for the upper and lower surfaces, trailing-edge spacing values $(Δz_u,~Δz_l)$, and a coefficient for the leading edge modifications at the nose.

source
AeroFuse.AircraftGeometry.camber_coordinatesFunction
camber_coordinates(wing :: WingMesh, n_span = wing.num_span, n_chord = wing.num_chord)

Generate the camber coordinates of a WingMesh with default spanwise $n_s$ and chordwise $n_c$ panel distributions from the mesh.

source
AeroFuse.AircraftGeometry.camber_thicknessFunction
camber_thickness(foil :: Foil, num :: Integer)

Compute the camber-thickness distribution of a Foil with cosine interpolation. Optionally specify the number of points for interpolation, default is 40.

source
AeroFuse.AircraftGeometry.camber_thicknessMethod
camber_thickness(wing :: Wing, num :: Integer)

Compute the camber-thickness distribution at each spanwise intersection of a Wing. A num must be specified to interpolate the internal Foil coordinates, which affects the accuracy of $(t/c)ₘₐₓ$ accordingly.

source
AeroFuse.AircraftGeometry.camber_thickness_to_CSTMethod
camber_thickness_to_CST(coords, num_dvs)

Convert camber-thickness coordinates to a specified number of Bernstein polynomial coefficients under a Class Shape Transformation by performing a least-squares solution.

source
AeroFuse.AircraftGeometry.chord_coordinatesFunction
chord_coordinates(wing :: WingMesh, n_span = wing.num_span, n_chord = wing.num_chord)

Generate the chord coordinates of a WingMesh with default spanwise $n_s$ and chordwise $n_c$ panel distributions from the mesh.

source
AeroFuse.AircraftGeometry.control_surfaceMethod
control_surface(foil :: Foil, δ, xc_hinge)
control_surface(foil :: Foil; angle, hinge)

Modify a Foil to mimic a control surface by specifying a deflection angle $δ$(in degrees, clockwise-positive convention) and a normalized hingex-coordinate∈ [0,1]in terms of the chord length. A constructor with named argumentsangle, hinge` is provided for convenience.

source
AeroFuse.AircraftGeometry.coordinatesFunction
wetted_area(fuse :: HyperEllipseFuselage, t)

Get the coordinates of a HyperEllipseFuselage given the parameter distribution $t$. Note that the distribution must have endpoints 0 and 1.

source
AeroFuse.AircraftGeometry.kulfan_CSTFunction
kulfan_CST(alpha_u, alpha_l,
           (Δz_u, Δz_l) = (0., 0.),
           (LE_u, LE_l) = (0., 0.),
           n            = 40)

Define a cosine-spaced foil with $2n$ points using the Class Shape Transformation method on a Bernstein polynomial basis for the upper and lower coordinates.

The foil is defined by arrays of coefficients $(α_u,~ α_l)$ for the upper and lower surfaces (not necessarily of the same lengths), trailing-edge displacement values $(Δz_u,~ Δz_l)$, and coefficients for leading edge modifications on the upper and lower surfaces at the nose.

source
AeroFuse.AircraftGeometry.maximum_thickness_to_chordFunction
maximum_thickness_to_chord(foil :: Foil, num :: Integer)

Compute the maximum thickness-to-chord ratio $(t/c)ₘₐₓ$ and its location $(x/c)$ of a Foil. Returned as the pair $(x/c, (t/c)ₘₐₓ)$.

A num must be specified to interpolate the Foil coordinates, which affects the accuracy of $(t/c)ₘₐₓ$ accordingly, default is 40.

source
AeroFuse.AircraftGeometry.maximum_thickness_to_chordMethod
maximum_thickness_to_chord(wing :: Wing, num :: Integer)

Compute the maximum thickness-to-chord ratios $(t/c)ₘₐₓ$ and their locations $(x/c)$ at each spanwise intersection of a Wing.

Returns an array of pairs $[(x/c),(t/c)ₘₐₓ]$, in which the first entry of each pair is the location (normalized to the local chord length at the spanwise intersection) and the corresponding maximum thickness-to-chord ratio at the intersection.

A num must be specified to interpolate the internal Foil coordinates, which affects the accuracy of $(t/c)ₘₐₓ$ accordingly.

source
AeroFuse.AircraftGeometry.mean_aerodynamic_centerFunction
mean_aerodynamic_center(wing :: Wing, 
    factor = 0.25; 
    symmetry = wing.symmetry, 
    flip = wing.flip
)

Compute the mean aerodynamic center of a Wing. By default, the factor is assumed to be at 25% from the leading edge, which can be adjusted. Similarly, options are provided to account for symmetry or to flip the location in the $x$-$z$ plane.

source
AeroFuse.AircraftGeometry.naca4Function
naca4(digits :: NTuple{4, <: Real}, n :: Integer = 40; sharp_TE :: Bool)
naca4(a, b, c, d, n :: Integer = 40; sharp_TE :: Bool)

Generate a Foil of a NACA 4-digit series profile with the specified digits, number of points (40 by default), and a named option to specify a sharp or blunt trailing edge.

Refer to the formula for the digits here: http://airfoiltools.com/airfoil/naca4digit

source
AeroFuse.AircraftGeometry.naca4_coordinatesMethod
naca4_coordinates(digits :: NTuple{4, <: Real}, n :: Integer, sharp_TE :: Bool)

Generate the coordinates of a NACA 4-digit series profile with a specified number of points, and a Boolean flag to specify a sharp or blunt trailing edge.

source
AeroFuse.AircraftGeometry.read_foilMethod
read_foil(
    path :: String; 
    header = true; 
    name = ""
)

Generate a Foil from a file consisting of 2D coordinates with named arguments to skip the header (first line of the file) or assign a name.

By default, the header is assumed to exist and should contain the airfoil name, which is assigned to the name of the Foil.

source
AeroFuse.AircraftGeometry.rotateMethod
rotate(foil :: Foil; 
    angle, 
    center = zeros(2)
)

Rotate the coordinates of a Foil about a 2-dimensional point (default is origin) by the angle $θ$ (in degrees).

source
AeroFuse.AircraftGeometry.surface_coordinatesFunction
surface_coordinates(wing :: WingMesh, n_span = wing.num_span, n_chord = wing.num_chord)

Generate the surface coordinates of a WingMesh with default spanwise $n_s$ and chordwise $n_c$ panel distributions from the mesh.

source
AeroFuse.AircraftGeometry.surface_panelsFunction
surface_panels(
    wing_mesh :: WingMesh, 
    n_s = wing_mesh.num_span, 
    n_c = length(first(foils(wing_mesh.surface))).x
)

Generate the surface panel distribution from a WingMesh with the default spanwise $n_s$ panel distribution from the mesh and the chordwise panel $n_c$ distribution from the number of points defining the root airfoil.

In case of strange results, provide a higher number of chordwise panels to represent the airfoils more accurately

source
AeroFuse.AircraftGeometry.taper_ratioMethod
taper_ratio(wing :: Wing)

Compute the taper ratio of a Wing, defined as the tip chord length divided by the root chord length, independent of the number of sections.

source
AeroFuse.AircraftGeometry.volumeMethod
volume(fuse :: HyperEllipseFuselage, ts)

Compute the volume of a HyperEllipseFuselage given the parameter distribution $t$. Note that the distribution must have endpoints 0 and 1.

source
AeroFuse.AircraftGeometry.wetted_area_ratioFunction
wetted_area_ratio(
    wing_mesh :: WingMesh, 
    n_s = wing_mesh.num_span, 
    n_c = length(first(foils(wing_mesh.surface))).x
)

Determine the wetted area ratio $S_{wet}/S$ of a WingMesh by calculating the ratio of the total area of the surface panels to the projected area of the Wing.

The wetted area ratio should be slightly above 2 for thin airfoils.

source
AeroFuse.PanelGeometry.make_panelsMethod
make_panels(foil :: Foil)
make_panels(foil :: Foil, n :: Integer)

Generate a vector of Panel2Ds from a Foil, additionally with cosine interpolation using (approximately) $n$ points if provided.

source
AeroFuse.PanelGeometry.wetted_areaFunction
wetted_area(
    wing_mesh :: WingMesh, 
    n_s = wing_mesh.num_span, 
    n_c = length(first(foils(wing_mesh.surface))).x
)

Determine the wetted area $S_{wet}$ of a WingMesh by calculating the total area of the surface panels.

source
AeroFuse.PanelGeometry.wetted_areaMethod
wetted_area(fuse :: HyperEllipseFuselage, t)

Compute the wetted area of a HyperEllipseFuselage given the parameter distribution $t$. Note that the distribution must have endpoints 0 and 1.

source

Panel Geometry API

AeroFuse.PanelGeometry.Panel3DType
Panel3D(p1, p2, p3, p4)

Four Cartesian coordinates p1, p2, p3, p4 representing corners of a panel in 3 dimensions. The following commutative diagram (math joke) depicts the order:

z → y
↓
x
        p1 —→— p4
        |       |
        ↓       ↓
        |       |
        p2 —→— p3
source
AeroFuse.PanelGeometry.transformMethod
transform(panel :: Panel3D, rotation, translation)

Perform an affine transformation on the coordinates of a Panel3D given a rotation matrix and translation vector.

source
AeroFuse.PanelGeometry.transform_normalMethod
transform_normal(panel :: Panel3D, h_l, g_l)

Transform the normal vector $n̂₀$ of a Panel3D about the hinge axis $ĥₗ$ by the control gain $gₗ$.

The transformation is the following: $n̂ₗ = gₗ ĥₗ × n̂₀$ (Flight Vehicle Aerodynamics, M. Drela, Eq. 6.36).

source
AeroFuse.PanelGeometry.wake_panelMethod
wake_panel(panels :: DenseArray{<: AbstractPanel3D}, bound, U)

Calculate required transformation from the global coordinate system to an to an AbstractPanel3D's local coordinate system.

source