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
AeroFuse.PanelGeometry.get_transformationFunction
get_transformation(panel :: AbstractPanel3D)

Generate the mapping to transform a point from global coordinates to an AbstractPanel3D's local coordinate system.

AeroFuse.PanelGeometry.make_panelsMethod
make_panels(xyzs)

Convert an array of coordinates corresponding to a wing, ordered from root to tip and leading-edge to trailing-edge, into panels.

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.

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).

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.

AeroFuse.DoubletSource.doublet_matrixMethod
doublet_matrix(panels_1, panels_2)

Create the matrix of doublet potential influence coefficients between pairs of panels₁ and panels₂.

AeroFuse.DoubletSource.influence_matrixMethod
influence_matrix(panels, wake_panel :: AbstractPanel2D)

Assemble the Aerodynamic Influence Coefficient matrix consisting of the doublet matrix, wake vector, Kutta condition given Panel2Ds and the wake panel.

AeroFuse.DoubletSource.solve_linearMethod
solve_linear(panels, u, wakes)

Solve the linear aerodynamic system given the array of Panel2Ds, a velocity $\vec U$, a vector of wake Panel2Ds, and an optional named bound for the length of the wake.

The system of equations $A[φ] = [\vec{U} ⋅ n̂] - B[σ]$ is solved, where $A$ is the doublet influence coefficient matrix, $φ$ is the vector of doublet strengths, $B$ is the source influence coefficient matrix, and $σ$ is the vector of source strengths.

AeroFuse.DoubletSource.solve_linearMethod
solve_linear(panels, u, sources, bound)

Solve the system of equations $[AIC][φ] = [\vec{U} ⋅ n̂] - B[σ]$ condition given the array of Panel2Ds, a velocity $\vec U$, a condition whether to disable source terms ($σ = 0$), and an optional named bound for the length of the wake.

AeroFuse.DoubletSource.source_strengthsMethod
source_strengths(panels, freestream)

Create the vector of source strengths for the Dirichlet boundary condition $σ = \vec U_{\infty} ⋅ n̂$ given Panel2Ds and a Uniform2D.

AeroFuse.DoubletSource.wake_vectorMethod
wake_vector(woke_panel :: AbstractPanel2D, panels)

Create the vector of doublet potential influence coefficients from the wake on the panels given the wake panel and the array of Panel2Ds.

AeroFuse.surface_velocitiesMethod
surface_velocities(panels, φs, u, sources :: Bool)

Compute the surface speeds and panel distances given the array of Panel2Ds, their associated doublet strengths $φ$s, the velocity $u$, and a condition whether to disable source terms ($σ = 0$).

AeroFuse.NonDimensional.aerodynamic_coefficientsMethod
aerodynamic_coefficients(force, moment, Ω, V, S, b, c, ρ)

Compute the relevant aerodynamic coefficients given a net force, moment and angular rates $Ω$ with reference speed $V$, area $S$, span $b$, chord $c$, and density $ρ$.

AeroFuse.NonDimensional.rate_coefficientMethod
rate_coefficient(Ω, V, L)

Compute the non-dimensional angular velocity coefficient given angular speed $Ω$, reference speed $V$, and reference length $L$.

AeroFuse.Beams.MaterialMethod
Material(E, J, σ_max, ρ)
Material(; elastic_modulus, shear_modulus, yield_stress, density)

Define a Material with positive, real elastic modulus $E$, shear modulus $G$, yield stress $σ_{\max}$, and density $ρ$.

The default assignments for the named variables are set to the properties of aluminium?

Arguments

  • elastic_modulus :: Real = 85e9: Elastic modulus.
  • shear_modulus :: Real = 25e9: Shear modulus.
  • yield_stress :: Real = 350e6: Yield stress.
  • density :: Real = 1.6e3: Density.
AeroFuse.Beams.TubeType
Tube(material :: Material, length, radius, thickness)

Define a hollow tube of fixed radius with a given material, length, and thickness.

AeroFuse.Beams.tube_stiffness_matrixFunction
tube_stiffness_matrix(E, G, A, Iyy, Izz, J, L)
tube_stiffness_matrix(E, G, A, Iyy, Izz, J, L, num)
tube_stiffness_matrix(tube :: Tube)
tube_stiffness_matrix(tube :: Vector{Tube})
tube_stiffness_matrix(x :: Matrix{Real})

Generate the stiffness matrix using the properties of a tube.

The required properties are the elastic modulus $E$, shear modulus $G$, area $A$, moments of inertia about the $y-$ and $z-$ axes $I_{yy}, I_{zz}$, polar moment of inertia $J$, length $L$. A composite stiffness matrix is generated with a specified number of elements.

AeroFuse.Laplace.FreestreamType
Freestream(α, β, Ω)
Freestream(α, β, Ω_x, Ω_y, Ω_z)
Freestream(U, Ω)
Freestream(; 
    alpha = 0., 
    beta  = 0., 
    omega = [0,0,0]
)

Define freestream conditions with angle of attack $α$ (degrees), sideslip angle $β$ (degrees), and (quasi-steady) rotation vector $Ω$ (in geometry axes).

Alternatively, provide the velocity vector $U$, which is normalized to determine the angles.

AeroFuse.velocityMethod
velocity(freestream :: Freestream)

Compute the velocity of a Freestream.

AeroFuse.VortexLattice.HorseshoeType
Horseshoe(panel :: Panel3D, normal, drift = zeros(3))

Generate a Horseshoe corresponding to a Panel3D, an associated normal vector, and a "drift velocity".

AeroFuse.VortexLattice.VortexRingMethod

Constructor for making a VortexRing with a Panel3D. The following convention is adopted:

    p1 —front leg→ p4
    |               |
left leg       right leg
    ↓               ↓
    p2 —back leg-→ p3
AeroFuse.VortexLattice.control_pointMethod
control_point(panel :: Panel3D)

Compute the control point of a Panel3D for horseshoes/vortex rings, which is the average of the 3-quarter point on each side of the trailing legs.

AeroFuse.bound_legMethod
bound_leg(panel :: Panel3D)

Compute the bound leg for a Panel3D, for horseshoes/vortex rings, which quarter point on each side of the trailing legs.

AeroFuse.form_factorMethod
form_factor(fuse :: HyperEllipseFuselage)

Compute the form factor of HyperEllipseFuselage for parasitic drag estimation based on the DATCOM formula (Raymer - Aircraft Design: A Conceptual Approach, 4th Ed., Eq. 12.31): $FF = 1 + 60 / f³ + f / 400$

AeroFuse.form_factorMethod
form_factor(wing :: Wing, M, num)
form_factor(wing_mesh :: WingMesh, M, num)

Compute the form factor $K_f$ of a Wing or WingMesh at Mach number $M$ for parasitic drag estimation based on Raymer's formula for wing, tail, strut, and pylon (Raymer - Aircraft Design: A Conceptual Approach, 4th Ed., Eq. 12.30): $FF = (1 + 0.6(t/c)ₘₐₓ / (x/c) + 100(t/c)ₘₐₓ^4) ⨯ 1.34M^0.18 × cos(Λ_m)^0.28$.

An integer num is required for interpolating Wing's airfoils' coordinates to determine their maximum thickness-to-chord ratios.

AeroFuse.lateral_stability_derivativesMethod
lateral_stability_derivatives(dvs, U, m, Iyy, Q, S, c)

Compute the stability derivatives for the lateral forces and moments with respect to the sideslip angle, roll and yaw rates: $[Y,L,N]_{β,p,r}$

The inputs are force and moment coefficient stability derivatives matrix dvs, the freestream speed U, the mass m and lateral moment sof inertia $I_{xx}, I_{zz}$, the dynamic pressure Q, reference area $S$ and span length $b$.

AeroFuse.longitudinal_stability_derivativesFunction
longitudinal_stability_derivatives(dvs, U₀, m, Iyy, Q, S, c)

Compute the stability derivatives for the longitudinal forces and moments with resepct to speed, angle of attack, and pitch rate: $[X,Z,M]_{u,α,q}$

The inputs are force and moment coefficient stability derivatives matrix dvs, the freestream speed U, the mass m and longitudinal moment of inertia $I_{yy}$, the dynamic pressure Q, reference area $S$ and chord length $c$.

AeroFuse.make_horseshoesMethod
make_horseshoes(wing :: WingMesh)

Generate an array of Horseshoes defined by the chord coordinates and normal vectors defined by the camber distribution of a WingMesh.

AeroFuse.make_vortex_ringsMethod
make_vortex_rings(wing :: WingMesh)

Generate an array of VortexRings defined by the camber coordinates and normal vectors of a WingMesh.

AeroFuse.parasitic_drag_coefficientFunction
parasitic_drag_coefficient(
    fuse :: HyperEllipseFuselage, 
    refs :: References, 
    x_tr :: Real, 
    ts = 0:0.01:1
)

Estimate the profile drag coefficient of a HyperEllipseFuselage using the wetted-area method based on Schlichting's skin-friction coefficient formula with given References, a specified transition location $xₜᵣ$ as a ratio of the fuselage length, and optionally the parametric distribution t ∈ [0,1] for the discretization of the nose, cabin and rear sections.

AeroFuse.parasitic_drag_coefficientMethod
parasitic_drag_coefficient(
    L, x_tr, 
    ρ, V, M, μ, S_ref, 
    S_wet, Kf, fM
)

Compute the parasitic drag coefficient $C_{D₀}$ with the following quantities. Uses a transition-based model from laminar to turbulent flow based on Schlichting's skin-friction coefficient formulas.

Arguments

  • L = Reference length (m)
  • x_tr = Transition location as ratio of reference length
  • ρ = Density (kg/m³)
  • V = Speed (m/s)
  • M = Mach number
  • S_ref = Reference area (m²)
  • μ = Dynamic viscosity (kg/(m ⋅ s))
  • S_wet = Wetted area (m²)
  • Kf = Form factor
  • fM = Mach number influence
AeroFuse.parasitic_drag_coefficientMethod
parasitic_drag_coefficient(
    wing :: Wing, 
    refs :: References,
    x_tr
)

Estimate the profile drag coefficient of a Wing using the wetted-area method based on Schlichting's skin-friction coefficient formula with given References, a specified transition location $xₜᵣ ∈ [0,1]$ as a ratio of the chord lengths.

AeroFuse.parasitic_drag_coefficientMethod
parasitic_drag_coefficient(
    wing :: WingMesh, 
    refs :: References
    x_tr :: Real, 
    u_es, 
)

Estimate the profile drag coefficient of a WingMesh using the local-friction and local-dissipation method based on Schlichting's skin-friction coefficient formula with given References, a specified transition ratio $xₜᵣ$, and edge velocities $\mathbf u_e$.

At present, the edge velocities would be computed using the vortex lattice method via VortexLatticeSystem. For this case, the panels corresponding to the camber distribution are used in the calculation.

AeroFuse.plot_panelsMethod
plot_panels(panels :: Array{Panel3D})

Get the coordinates of an array of Panel3D for plotting.

AeroFuse.plot_planformMethod
plot_planform(wing :: AbstractWing)

Get the planform coordinates of an AbstractWing for plotting.

AeroFuse.plot_surfaceMethod
plot_surface(wing :: AbstractWing)
plot_surface(wing :: AbstractWing, 
             span_num :: Vector{Integer}, 
             chord_num :: Integer)

Get the surface coordinates of an AbstractWing for plotting, optionally with specified spanwise and chordwise discretizations.

AeroFuse.print_infoFunction
print_info(wing :: AbstractWing)

Print the relevant geometric characteristics of an AbstractWing.

AeroFuse.solve_caseMethod
solve_case(
    components :: Vector{Horseshoe}, 
    fs :: Freestream, 
    refs :: References;
    name = :aircraft :: Symbol, 
    compressible = false :: Boolean,
    print = false :: Boolean,
    print_components = false :: Boolean
)

Perform a vortex lattice analysis given a vector of Horseshoes, a Freestream condition, and References values.

AeroFuse.solve_caseMethod
solve_case(foil :: Foil, uniform :: Uniform2D; sources = false, wake_length = 1e3, num_panels :: Int64 = 60)

Evaluate a doublet-source case given a Foil with a Uniform2D, with optional named arguments to specify whether the source terms are non-zero, the length of the wake, and the number of panels for the analysis.

AeroFuse.spanwise_loadingMethod
spanwise_loading(panels :: Matrix{Panel3D}, CFs, S)
spanwise_loading(wing :: WingMesh, CFs, S)

Obtain the spanwise aerodynamic loads (CDi, CY, CL) for a given matrix of Panel3Ds, surface coefficients CFs over the panels, and reference area $S$.

If a WingMesh is provided instead of the matrix, then it will compute the chordwise panel distribution over the surface.

AeroFuse.spanwise_loadingMethod
spanwise_loading(wing_mesh :: WingMesh, ref :: References, CFs, Γs)

Obtain the spanwise aerodynamic loads (CDi, CY, CL, CL_norm) for a given WingMesh, reference values ref, surface coefficients CFs, and circulations $Γ$.

The chordwise normalized lift coefficient loading CL_norm is calculated via the formula $C_L = 2Γ / ρVc$ based on the reference speed $V$ and chord $c$ from ref.

Note that the surface coefficients CFs should be computed in wind axes to match CL_norm.

AeroFuse.MathTools.weightedMethod
weighted(x1, x2, μ)

Compute the weighted value between two values $x_1$ and $x_2$ with weight $\mu \in [0,1]$.

AeroFuse.VortexLattice.HorseshoeType
Horseshoe(r1, r2, rc, normal, chord)

Define a horseshoe vortex with a start and endpoints $r₁, r₂$ for the bound leg, a collocation point $r$, a normal vector $n̂$, and a finite core size.

The finite core setup is not implemented for now.

AeroFuse.VortexLattice.ReferencesType
References(V, ρ, μ, S, b, c, r)
References(; 
    speed, density, viscosity,
    sound_speed, area, span, 
    chord, location
)

Define reference values with speed $V$, density $ρ$, dynamic viscosity $μ$, area $S$, span $b$, chord $c$, location $r$ for a vortex lattice analysis. A constructor with named arguments is provided for convenience:

Arguments

  • speed :: Real = 1.: Speed (m/s)
  • density :: Real = 1.225: Density (m)
  • viscosity :: Real = 1.5e-5: Dynamic viscosity (kg/(m ⋅ s))
  • sound_speed :: Real = 330.: Speed of sound (m/s)
  • area :: Real = 1.: Area (m²)
  • span :: Real = 1.: Span length (m)
  • chord :: Real = 1.: Chord length (m)
  • location :: Vector{Real} = [0,0,0]: Position (m)
AeroFuse.VortexLattice.VortexLatticeSystemType
VortexLatticeSystem(
    aircraft, 
    fs :: Freestream, 
    refs :: References, 
    compressible = false, 
)

Construct a VortexLatticeSystem for analyzing inviscid aerodynamics of an aircraft (must be a ComponentArray of Horseshoes or VortexRings) with Freestream conditions and References for non-dimensionalization. Options are provided for compressibility corrections via the Prandtl-Glauert transformation (false by default) and axis system for computing velocities and forces (Geometry by default).

AeroFuse.VortexLattice.VortexLatticeSystemType
VortexLatticeSystem

A system consisting of the relevant variables for a vortex lattice analysis for post-processing.

Arguments

The accessible fields are:

  • vortices: The array of vortices, presently of AbstractVortex types.
  • circulations: The circulation strengths of the vortices obtained by solving the linear system.
  • influence_matrix: The influence matrix of the linear system.
  • boundary_vector: The boundary condition corresponding to the right-hand-side of the linear system.
  • freestream :: Freestream: The freestream conditions.
  • reference :: References: The reference values.
AeroFuse.VortexLattice.VortexRingType
VortexRing(r1, r2, r3, r4, r_c, n̂, ε)

A vortex ring consisting of four points $r_i, i = 1,…,4$, a collocation point $r_c$, a normal vector $n̂$, and a core size $ε$. The following convention is adopted:

    r1 —front leg→ r4
    |               |
left leg       right leg
    ↓               ↓
    r2 —back leg-→ r3
AeroFuse.VortexLattice.bound_velocityMethod
bound_velocity(r, hs :: Horseshoe, Γ, u_hat)

Compute the induced velocity at a point $r$ from the bound leg with constant strength $Γ$ of a given Horseshoe.

AeroFuse.VortexLattice.boundary_conditionMethod
boundary_condition(vortices, U, Ω)

Assemble the boundary condition vector given an array of AbstractVortex, the freestream velocity $U$, and a quasi-steady rotation vector $Ω$.

AeroFuse.VortexLattice.center_of_pressureMethod
center_of_pressure(system :: VortexLatticeSystem)

Determine the center of pressure $x_{cp}$ of the VortexLatticeSystem.

This is computed based on the nearfield lift $C_L$ and moment $Cₘ$ coefficients, and the reference location $xᵣ$ and chord length $cᵣ$ from References: $x_{cp} = xᵣ -cᵣ(Cₘ / C_L)$

AeroFuse.VortexLattice.farfieldMethod
farfield(system :: VortexLatticeSystem)

Compute the total farfield force coefficients of the VortexLatticeSystem. These are in wind axes by definition.

AeroFuse.VortexLattice.farfield_coefficientsMethod
farfield_coefficients(system :: VortexLatticeSystem)

Compute the total farfield force coefficients for all components of the VortexLatticeSystem. These are in wind axes by definition.

AeroFuse.VortexLattice.freestream_derivativesMethod
freestream_derivatives(
    system :: VortexLatticeSystem;
    axes :: AbstractAxisSystem = Stability(),
    name = :aircraft,
    print = false,
    print_components = false,
    farfield = false
)

Compute the force and moment coefficients of the components of a VortexLatticeSystem and their derivatives with respect to freestream values: Mach $M$ (if compressible), angles of attack $α$ and sideslip $β$, and non-dimensionalized rotation rates $p̄, q̄, r̄$ (in stability axes).

The axes of the force and moment coefficients can be changed by passing any AbstractAxisSystem (such as Body(), Geometry(), Wind(), Stability()) to the named axes argument. The nearfield force and moment coefficients are reported in stability axes by default. Note that the derivatives with respect to the rotation rates will still refer to the rotation vector in stability axes $(p̄, q̄, r̄)$.

Optional printing arguments are provided for the components and the entire system, along with the corresponding farfield coefficients if needed.

AeroFuse.VortexLattice.influence_matrixMethod
influence_matrix(vortices)
influence_matrix(vortices, u_hat)

Assemble the Aerodynamic Influence Coefficient (AIC) matrix given an array of AbstractVortex and the freestream direction $û$.

AeroFuse.VortexLattice.nearfieldMethod
nearfield(system :: VortexLatticeSystem)

Compute the total nearfield force and moment coefficients for all components of the VortexLatticeSystem. These are in wind axes by default.

AeroFuse.VortexLattice.nearfield_coefficientsMethod
nearfield_coefficients(system :: VortexLatticeSystem)

Compute the nearfield force and moment coefficients for all components of the VortexLatticeSystem. These are in wind axes by default.

AeroFuse.VortexLattice.print_coefficientsFunction
print_coefficients(
    system :: VortexLatticeSystem, 
    name = :aircraft;
    components = false
)

Print a pretty table of the total nearfield and farfield coefficients of a VortexLatticeSystem with an optional name.

A named Boolean argument components is provided to also enable the printing of any possible components.

AeroFuse.VortexLattice.print_derivativesFunction
print_derivatives(
    nf_coeffs, ff_coeffs, name = "";
    farfield = false,
)

Print a pretty table of the aerodynamic coefficients and derivatives with an optional name, and a named argument for enabling printing of farfield coefficients.

AeroFuse.VortexLattice.surface_dynamicsMethod
surface_dynamics(
    system :: VortexLatticeSystem; 
    axes :: AbstractAxisSystem = Geometry()
)

Compute the forces and moments for all components of the VortexLatticeSystem in a specified reference axis system as a named argument.

The reference axis system is set to the geometry axes defined in the construction of the VortexLatticeSystem by default.

AeroFuse.VortexLattice.surface_forcesMethod
surface_forces(
    system :: VortexLatticeSystem; 
    axes :: AbstractAxisSystem = Geometry()
)

Compute the forces for all components of the VortexLatticeSystem in a specified reference axis system as a named argument.

The reference axis system is set to the geometry axes defined in the construction of the VortexLatticeSystem by default.

AeroFuse.VortexLattice.surface_momentsMethod
surface_moments(
    system :: VortexLatticeSystem; 
    axes :: AbstractAxisSystem = Geometry()
)

Compute the moments for all components of the VortexLatticeSystem in a specified reference axis system as a named argument.

The reference axis system is set to the geometry axes defined in the construction of the VortexLatticeSystem by default.

AeroFuse.VortexLattice.trailing_velocityMethod
trailing_velocity(r, hs :: Horseshoes, Γ, u_hat)

Compute the induced velocity at a point $r$ from the semi-infinite trailing legs with constant strength $Γ$ of a given Horseshoe hs.

AeroFuse.VortexLattice.transformMethod
transform(hs :: Horseshoe, T :: LinearMap)

Generate a new Horseshoe with the points and normal vectors transformed by the LinearMap $T$.

AeroFuse.VortexLattice.transformMethod
transform(ring :: VortexRing, T :: LinearMap)

Generate a new VortexRing with the points and normal vectors transformed by the LinearMap $T$.

AeroFuse.solve_linearMethod
solve_linear(horseshoes, normals, U, Ω)

Evaluate and return the vortex strengths $Γ$s given an array of Horseshoes, their associated normal vectors, the velocity vector $U$, and the quasi-steady rotation vector $Ω$.

AeroFuse.surface_coefficientsMethod
surface_coefficients(
    system :: VortexLatticeSystem; 
    axes :: AbstractAxisSystem = Geometry()
)

Compute the force and moment coefficients of the surfaces over all components in a given VortexLatticeSystem, in a specified reference axis system as a named argument.

The reference axis system is set to the geometry axes defined in the construction of the VortexLatticeSystem by default.

AeroFuse.surface_velocitiesMethod
surface_velocities(
    system :: VortexLatticeSystem; 
    axes :: AbstractAxisSystem = Geometry()
)

Compute the induced velocities for all components of the VortexLatticeSystem in a specified reference axis system as a named argument.

The reference axis system is set to the geometry axes defined in the construction of the VortexLatticeSystem by default.

AeroFuse.velocityFunction
velocity(r, hs :: Horseshoe, Γ, u_hat = [1.,0.,0.])

Compute the induced velocity at a point $r$ of a given Horseshoe with a bound leg of constant strength $Γ$ and semi-infinite trailing legs pointing in a given direction $û$, by default û = x̂.

AeroFuse.velocityFunction
velocity(r, ring :: VortexRing, Γ)

Computes the velocity at a point $r$ induced by a VortexRing with constant strength $Γ$.

AeroFuse.velocityMethod
velocity(freestream :: Freestream, ::Geometry)

Compute the velocity of Freestream in the geometry axis system.

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.

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)
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)
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.

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)
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 $θ$.

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.

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.

AeroFuse.AircraftGeometry.camber_lineFunction
camber_line(foil :: Foil, n :: Integer = 40)

Get the camber line of a Foil. Optionally specify the number of points for linear interpolation, default is 40.

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.

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.

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.

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.

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.

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.

AeroFuse.AircraftGeometry.coordinates_to_CSTMethod
coordinates_to_CST(coords, num_dvs)

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

AeroFuse.AircraftGeometry.cosine_interpolationFunction
cosine_interpolation(foil :: Foil, n :: Integer = 40)

Interpolate a Foil profile's coordinates to a cosine by projecting the x-coordinates of a circle onto the geometry with $2n$ points.

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.

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.

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.

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.

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

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.

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.

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).

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.

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

AeroFuse.AircraftGeometry.sweepsFunction
sweeps(wing :: Wing, w = 0.)

Obtain the sweep angles (in radians) at the corresponding normalized chord length ratio $w ∈ [0,1]$.

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.

AeroFuse.AircraftGeometry.thickness_lineFunction
thickness_line(foil :: Foil, n :: Integer = 40)

Get the thickness line of a Foil. Optionally specify the number of points for linear interpolation, default is 40.

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.

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.

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.

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.

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.