Objectives

Here we will show you how to perform an aerodynamic stability analysis of a conventional aircraft. Here, we'll attempt to replicate the design of a Boeing 777. It won't be a realistic replication, but the overall geometry will match to a certain extent.

Source: boeingboeing2, DeviantArt

Note

Refer to the Aircraft Aerodynamic Analysis tutorial before studying this tutorial.

Recipe

  1. Define the geometries of a fuselage, wing, horizontal tail and vertical tail.
  2. Perform an aerodynamic analysis of this aircraft configuration at given freestream conditions and reference values.
  3. Evaluate the derivatives of the aerodynamic coefficients with respect to the freestream conditions.
  4. Evaluate the quantities of interest for aerodynamic stability.

Let's import the relevant packages.

using AeroFuse          # Main package
using Plots             # Plotting library
gr(                     # Plotting backend
    size = (800,600),   # Size
    dpi = 300,          # Resolution
    palette = :Dark2_8  # Color scheme
)
using LaTeXStrings      # For LaTeX printing in plots

Aircraft Geometry

First, let's define the fuselage. Here we'll define it by combining a cylindrical definition for the cabin with hyperelliptical curves for the nose and rear.

# Fuselage definition
fuse = HyperEllipseFuselage(
    radius = 3.04,          # Radius, m
    length = 63.7,          # Length, m
    x_a    = 0.15,          # Start of cabin, ratio of length
    x_b    = 0.7,           # End of cabin, ratio of length
    c_nose = 2.0,           # Curvature of nose
    c_rear = 1.2,           # Curvature of rear
    d_nose = -0.5,          # "Droop" or "rise" of nose, m
    d_rear = 1.0,           # "Droop" or "rise" of rear, m
    position = [0.,0.,0.]   # Set nose at origin, m
)

# Compute geometric properties
ts = 0:0.1:1                # Distribution of sections for nose, cabin and rear
S_f = wetted_area(fuse, ts) # Surface area, m²
V_f = volume(fuse, ts)      # Volume, m³

# Get coordinates of rear end
fuse_end = fuse.affine.translation + [ fuse.length, 0., 0. ]
3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
 63.7
  0.0
  0.0

Now, let's define the lifting surfaces. We'll download a supercritical airfoil for the wing section; note that this is not the same one as used in the Boeing 777-200LR. We'll also define a two-section wing.

# Define one airfoil
foil_w = read_foil(download("http://airfoiltools.com/airfoil/seligdatfile?airfoil=rae2822-il"))

# Define vector of airfoils
foils = [ foil_w, foil_w, naca4((0,0,1,2)) ]

# Wing
wing = Wing(
    foils       = foils,                        # Airfoils
    chords      = [14.0, 9.73, 1.43561],        # Chord lengths
    spans       = [14.0, 46.9] / 2,             # Span lengths
    dihedrals   = fill(6, 2),                   # Dihedral angles (deg)
    sweeps      = fill(35.6, 2),                # Sweep angles (deg )
    w_sweep     = 0.,                           # Leading-edge sweep
    position    = [19.51, 0., -2.5],            # Position
    symmetry    = true                          # Symmetry
)

b_w = span(wing)
S_w = projected_area(wing)
c_w = mean_aerodynamic_chord(wing)
x_w, y_w, z_w = mac_w = mean_aerodynamic_center(wing)
3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
 30.683074339989055
  0.0
 -2.5

For reference, let's plot what we have so far.

p1 = plot(
    xaxis = L"x", yaxis = L"y", zaxis = L"z",
    aspect_ratio = 1,
    zlim = (-0.5, 0.5) .* span(wing),
    camera = (30,30)
)

plot!(fuse, label = "Fuselage", alpha = 0.6)
plot!(wing, label = "Wing")

Stabilizers

Now, let's add the stabilizers. First, the horizontal tail.

htail = WingSection(
    area        = 101,  # Area (m²)
    aspect      = 4.2,  # Aspect ratio
    taper       = 0.4,  # Taper ratio
    dihedral    = 7.,   # Dihedral angle (deg)
    sweep       = 35.,  # Sweep angle (deg)
    w_sweep     = 0.,   # Leading-edge sweep
    root_foil   = naca4(0,0,1,2),
    symmetry    = true,

    # Orientation
    angle       = -2,           # Incidence angle (deg)
    axis        = [0., 1., 0.], # Axis of rotation, y-axis
    position    = [ fuse_end.x - 8., 0., 0.],
)

b_h = span(htail)
S_h = projected_area(htail)
c_h = mean_aerodynamic_chord(htail)
x_h, y_h, z_h = mac_h = mean_aerodynamic_center(htail)
3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
 60.08867631466764
  0.0
  0.1532559539584015

Now the vertical tail.

vtail = WingSection(
    area        = 56.1, # Area (m²)
    aspect      = 1.5,  # Aspect ratio
    taper       = 0.4,  # Taper ratio
    sweep       = 44.4, # Sweep angle (deg)
    w_sweep     = 0.,   # Leading-edge sweep
    root_foil   = naca4(0,0,0,9),

    # Orientation
    angle       = 90.,       # To make it vertical
    axis        = [1, 0, 0], # Axis of rotation, x-axis
    position    = htail.affine.translation - [2.,0.,0.]
) # Not a symmetric surface

b_v = span(vtail)
S_v = projected_area(vtail)
c_v = mean_aerodynamic_chord(vtail)
x_v, y_v, z_v = mac_v = mean_aerodynamic_center(vtail)
3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
 59.17243218476657
  2.407305072320327e-16
  3.9314275332224553

Let's mesh and plot the lifting surfaces.

wing_mesh = WingMesh(wing, [8,16], 10,
    span_spacing = fill(Uniform(), 4) # Number of spacings = number of spanwise stations
    # (including symmetry)
)

htail_mesh = WingMesh(htail, [10], 8)
vtail_mesh = WingMesh(vtail, [8], 6)

# Plot meshes
plt = plot(
    xaxis = L"x", yaxis = L"y", zaxis = L"z",
    aspect_ratio = 1,
    zlim = (-0.5, 0.5) .* span(wing_mesh),
    camera = (30, 45),
)
plot!(fuse, label = "Fuselage", alpha = 0.6)
plot!(plt, wing_mesh, label = "Wing")
plot!(plt, htail_mesh, label = "Horizontal Tail")
plot!(plt, vtail_mesh, label = "Vertical Tail")

Aerodynamic Analysis

Now, let's generate the horseshoe system for the aircraft.

aircraft = ComponentVector(
    wing  = make_horseshoes(wing_mesh),
    htail = make_horseshoes(htail_mesh),
    vtail = make_horseshoes(vtail_mesh)
)
ComponentVector{Horseshoe{Float64}}(wing = Horseshoe{Float64}[Horseshoe{Float64}([41.31884181869616, -30.45, 0.7004239638398483], [39.226616003049074, -27.51875, 0.39233717421733394], [40.29663791737372, -28.984375, 0.5463805690285913], [0.0005456281159104126, 0.029861367098509173, 0.28033310122543664], 0.0) Horseshoe{Float64}([39.226616003049074, -27.51875, 0.39233717421733394], [37.13439018740198, -24.5875, 0.08425038459481993], [38.217098237408095, -26.053125, 0.23829377940607704], [0.0016368843477312378, 0.04627635308596634, 0.42907804209054934], 0.0) … Horseshoe{Float64}([37.13439018740198, 24.5875, 0.08425038459481993], [39.226616003049074, 27.51875, 0.39233717421733394], [38.217098237408095, 26.053125, 0.23829377940607704], [0.0016368843477312378, -0.04627635308596634, 0.42907804209054934], 0.0) Horseshoe{Float64}([39.226616003049074, 27.51875, 0.39233717421733394], [41.31884181869616, 30.45, 0.7004239638398483], [40.29663791737372, 28.984375, 0.5463805690285913], [0.0005456281159104126, -0.029861367098509173, 0.28033310122543664], 0.0); Horseshoe{Float64}([41.37067989644471, -30.45, 0.7004239638398485], [39.31589158416699, -27.51875, 0.39233717421733416], [40.412672379669004, -28.984375, 0.5463805690285911], [0.00019510184898996474, 0.08570486747568318, 0.8135583465329406], 0.0) Horseshoe{Float64}([39.31589158416699, -27.51875, 0.39233717421733416], [37.261103271889276, -24.5875, 0.08425038459481993], [38.39470066708566, -26.053125, 0.23829377940607693], [0.0005853055469684509, 0.13137909281537086, 1.2452329779493736], 0.0) … Horseshoe{Float64}([37.261103271889276, 24.5875, 0.08425038459481993], [39.31589158416699, 27.51875, 0.39233717421733416], [38.39470066708566, 26.053125, 0.23829377940607693], [0.0005853055469684509, -0.13137909281537086, 1.2452329779493736], 0.0) Horseshoe{Float64}([39.31589158416699, 27.51875, 0.39233717421733416], [41.37067989644471, 30.45, 0.7004239638398485], [40.412672379669004, 28.984375, 0.5463805690285911], [0.00019510184898996474, -0.08570486747568318, 0.8135583465329406], 0.0); … ; Horseshoe{Float64}([42.63406946277126, -30.45, 0.7004239638398482], [41.491702158378985, -27.51875, 0.39233717421733383], [42.13227244993827, -28.984375, 0.546380569028591], [0.02071501380863916, 0.09146069787357614, 0.8135583465329401], 0.0) Horseshoe{Float64}([41.491702158378985, -27.51875, 0.39233717421733383], [40.349334853986704, -24.5875, 0.08425038459481993], [41.02672174524036, -26.053125, 0.23829377940607688], [0.06214504142591892, 0.15166120081756396, 1.245232977949374], 0.0) … Horseshoe{Float64}([40.349334853986704, 24.5875, 0.08425038459481993], [41.491702158378985, 27.51875, 0.39233717421733383], [41.02672174524036, 26.053125, 0.23829377940607688], [0.06214504142591892, -0.15166120081756396, 1.245232977949374], 0.0) Horseshoe{Float64}([41.491702158378985, 27.51875, 0.39233717421733383], [42.63406946277126, 30.45, 0.7004239638398482], [42.13227244993827, 28.984375, 0.546380569028591], [0.02071501380863916, -0.09146069787357614, 0.8135583465329401], 0.0); Horseshoe{Float64}([42.7193199413754, -30.45, 0.7004239638398482], [41.6385206043654, -27.51875, 0.39233717421733383], [42.202829279371514, -28.984375, 0.546380569028591], [0.010666838939015766, 0.03320013990913556, 0.2803331012254997], 0.0) Horseshoe{Float64}([41.6385206043654, -27.51875, 0.39233717421733383], [40.5577212673554, -24.5875, 0.08425038459481993], [41.13471607804297, -26.053125, 0.23829377940607688], [0.03200051681704619, 0.05656144082582942, 0.4290780420905902], 0.0) … Horseshoe{Float64}([40.5577212673554, 24.5875, 0.08425038459481993], [41.6385206043654, 27.51875, 0.39233717421733383], [41.13471607804297, 26.053125, 0.23829377940607688], [0.03200051681704619, -0.05656144082582942, 0.4290780420905902], 0.0) Horseshoe{Float64}([41.6385206043654, 27.51875, 0.39233717421733383], [42.7193199413754, 30.45, 0.7004239638398482], [42.202829279371514, 28.984375, 0.546380569028591], [0.010666838939015766, -0.03320013990913556, 0.2803331012254997], 0.0)], htail = Horseshoe{Float64}[Horseshoe{Float64}([62.88890365446351, -10.298058069364341, 1.5162553193557646], [62.54031413176851, -9.794035232054846, 1.44215835982815], [62.76985874207378, -10.046046650709593, 1.4811362068318652], [-0.003889780601881568, 0.013685154467470186, 0.11138874253045755], 0.0) Horseshoe{Float64}([62.54031413176851, -9.794035232054846, 1.44215835982815], [61.528667934930915, -8.331303987175813, 1.2271206079077435], [62.09737450683494, -9.06266960961533, 1.3368354231503647], [-0.012848276001000303, 0.0452032285906627, 0.36792648581346676], 0.0) … Horseshoe{Float64}([61.528667934930915, 8.331303987175813, 1.2271206079077435], [62.54031413176851, 9.794035232054846, 1.44215835982815], [62.09737450683494, 9.06266960961533, 1.3368354231503647], [-0.012848276001000303, -0.0452032285906627, 0.36792648581346676], 0.0) Horseshoe{Float64}([62.54031413176851, 9.794035232054846, 1.44215835982815], [62.88890365446351, 10.298058069364341, 1.5162553193557646], [62.76985874207378, 10.046046650709593, 1.4811362068318652], [-0.003889780601881568, -0.013685154467470186, 0.11138874253045755], 0.0); Horseshoe{Float64}([63.04472763081975, -10.298058069364341, 1.5216968125153882], [62.707577960496174, -9.794035232054846, 1.4479993414354557], [63.03349105386857, -10.046046650709593, 1.4903424500226277], [-0.011077157969957813, 0.038972022690971254, 0.317208301282571], 0.0) Horseshoe{Float64}([62.707577960496174, -9.794035232054846, 1.4479993414354557], [61.7291315083166, -8.331303987175813, 1.2341209501454369], [62.39743171606405, -9.06266960961533, 1.3473136517881916], [-0.03658879445176172, 0.12872790398734174, 1.0477659852371668], 0.0) … Horseshoe{Float64}([61.7291315083166, 8.331303987175813, 1.2341209501454369], [62.707577960496174, 9.794035232054846, 1.4479993414354557], [62.39743171606405, 9.06266960961533, 1.3473136517881916], [-0.03658879445176172, -0.12872790398734174, 1.0477659852371668], 0.0) Horseshoe{Float64}([62.707577960496174, 9.794035232054846, 1.4479993414354557], [63.04472763081975, 10.298058069364341, 1.5216968125153882], [63.03349105386857, 10.046046650709593, 1.4903424500226277], [-0.011077157969957813, -0.038972022690971254, 0.317208301282571], 0.0); … ; Horseshoe{Float64}([65.32850451577563, -10.298058069364341, 1.6014480586855144], [65.15901884058357, -9.794035232054846, 1.533605543331635], [65.40109993639021, -10.046046650709593, 1.5730211740557802], [-0.011077157969957705, 0.03897202269097152, 0.31720830128257105], 0.0) Horseshoe{Float64}([65.15901884058357, -9.794035232054846, 1.533605543331635], [64.66715225377132, -8.331303987175813, 1.336718895359434], [65.0921625288351, -9.06266960961533, 1.4414157253432798], [-0.036588794451761386, 0.12872790398734185, 1.0477659852371666], 0.0) … Horseshoe{Float64}([64.66715225377132, 8.331303987175813, 1.336718895359434], [65.15901884058357, 9.794035232054846, 1.533605543331635], [65.0921625288351, 9.06266960961533, 1.4414157253432798], [-0.036588794451761386, -0.12872790398734185, 1.0477659852371666], 0.0) Horseshoe{Float64}([65.15901884058357, 9.794035232054846, 1.533605543331635], [65.32850451577563, 10.298058069364341, 1.6014480586855144], [65.40109993639021, 10.046046650709593, 1.5730211740557802], [-0.011077157969957705, -0.03897202269097152, 0.31720830128257105], 0.0); Horseshoe{Float64}([65.58280216782836, -10.298058069364341, 1.6103283283751404], [65.43198581212043, -9.794035232054846, 1.543137760023534], [65.56264383893216, -10.046046650709593, 1.5786624114392451], [-0.003889780601881672, 0.013685154467470526, 0.11138874253046113], 0.0) Horseshoe{Float64}([65.43198581212043, -9.794035232054846, 1.543137760023534], [64.99429970069268, -8.331303987175813, 1.3481431359431886], [65.27602622989178, -9.06266960961533, 1.4478363872657793], [-0.01284827600100058, 0.045203228590662406, 0.3679264858134668], 0.0) … Horseshoe{Float64}([64.99429970069268, 8.331303987175813, 1.3481431359431886], [65.43198581212043, 9.794035232054846, 1.543137760023534], [65.27602622989178, 9.06266960961533, 1.4478363872657793], [-0.01284827600100058, -0.045203228590662406, 0.3679264858134668], 0.0) Horseshoe{Float64}([65.43198581212043, 9.794035232054846, 1.543137760023534], [65.58280216782836, 10.298058069364341, 1.6103283283751404], [65.56264383893216, 10.046046650709593, 1.5786624114392451], [-0.003889780601881672, -0.013685154467470526, 0.11138874253046113], 0.0)], vtail = Horseshoe{Float64}[Horseshoe{Float64}([53.84630872675395, 0.0, 0.0], [54.18486988563193, 2.137860520751363e-17, 0.34913911868137404], [54.304865633097386, 1.0689302603756815e-17, 0.17456955934068702], [0.0, -0.40399072732320807, 2.47372975550789e-17], 0.0) Horseshoe{Float64}([54.18486988563193, 2.137860520751363e-17, 0.34913911868137404], [55.14901049489124, 8.225972198474937e-17, 1.3434032088602492], [54.9433606669902, 5.1819163596131497e-17, 0.8462711637708116], [0.0, -1.0993398152055094, 6.731514929333348e-17], 0.0) … Horseshoe{Float64}([61.439011047620205, 4.794447948899935e-16, 7.829927701992145], [62.40315165687952, 5.403259116672293e-16, 8.824191792171021], [62.05432531043232, 5.098853532786114e-16, 8.327059747081583], [0.0, -0.5299187314164414, 3.244816391186853e-17], 0.0) Horseshoe{Float64}([62.40315165687952, 5.403259116672293e-16, 8.824191792171021], [62.741712815757495, 5.617045168747429e-16, 9.173330910852394], [62.692820344325135, 5.510152142709861e-16, 8.998761351511707], [0.0, -0.16812879171659315, 1.0294919331011913e-17], 0.0); Horseshoe{Float64}([54.6849577820983, 0.0, 0.0], [55.00436743353075, 2.137860520751363e-17, 0.34913911868137404], [55.63498023034437, 1.0689302603756815e-17, 0.17456955934068702], [0.0, -1.1037231928337101, 6.7583553762425e-17], 0.0) Horseshoe{Float64}([55.00436743353075, 2.137860520751363e-17, 0.34913911868137404], [55.91396916385321, 8.225972198474937e-17, 1.3434032088602492], [56.21436308536693, 5.1819163596131497e-17, 0.8462711637708116], [0.0, -3.003452229924818, 1.8390840798847047e-16], 0.0) … Horseshoe{Float64}([61.84816105614035, 4.794447948899935e-16, 7.829927701992145], [62.75776278646279, 5.403259116672293e-16, 8.824191792171021], [62.66699118485762, 5.098853532786114e-16, 8.327059747081583], [0.0, -1.4477648981121707, 8.865003241954816e-17], 0.0) Horseshoe{Float64}([62.75776278646279, 5.403259116672293e-16, 8.824191792171021], [63.07717243789525, 5.617045168747429e-16, 9.173330910852394], [63.24637403988018, 5.510152142709861e-16, 8.998761351511707], [0.0, -0.45933640118489183, 2.8126242672147173e-17], 0.0); … ; Horseshoe{Float64}([60.65210209711992, 0.0, 0.0], [60.835245204102755, 2.137860520751363e-17, 0.34913911868137404], [61.53399127314118, 1.0689302603756815e-17, 0.17456955934068702], [0.0, -1.1037231928337075, 6.758355376242484e-17], 0.0) Horseshoe{Float64}([60.835245204102755, 2.137860520751363e-17, 0.34913911868137404], [61.3567926472095, 8.225972198474937e-17, 1.3434032088602492], [61.85121371233109, 5.1819163596131497e-17, 0.8462711637708116], [0.0, -3.0034522299248185, 1.8390840798847047e-16], 0.0) … Horseshoe{Float64}([64.75933961381432, 4.794447948899935e-16, 7.829927701992145], [65.28088705692106, 5.403259116672293e-16, 8.824191792171021], [65.38414259892375, 5.098853532786114e-16, 8.327059747081583], [0.0, -1.4477648981121634, 8.865003241954773e-17], 0.0) Horseshoe{Float64}([65.28088705692106, 5.403259116672293e-16, 8.824191792171021], [65.46403016390389, 5.617045168747429e-16, 9.173330910852394], [65.70136503811365, 5.510152142709861e-16, 8.998761351511707], [0.0, -0.4593364011848968, 2.812624267214748e-17], 0.0); Horseshoe{Float64}([61.99757944912139, 0.0, 0.0], [62.14999704659526, 2.137860520751363e-17, 0.34913911868137404], [62.36306457476276, 1.0689302603756815e-17, 0.17456955934068702], [0.0, -0.40399072732321056, 2.4737297555079053e-17], 0.0) Horseshoe{Float64}([62.14999704659526, 2.137860520751363e-17, 0.34913911868137404], [62.58404564147047, 8.225972198474937e-17, 1.3434032088602492], [62.64344182076147, 5.1819163596131497e-17, 0.8462711637708116], [0.0, -1.0993398152055023, 6.731514929333304e-17], 0.0) … Horseshoe{Float64}([65.4157549123554, 4.794447948899935e-16, 7.829927701992145], [65.84980350723059, 5.403259116672293e-16, 8.824191792171021], [65.76602316797545, 5.098853532786114e-16, 8.327059747081583], [0.0, -0.5299187314164555, 3.24481639118694e-17], 0.0) Horseshoe{Float64}([65.84980350723059, 5.403259116672293e-16, 8.824191792171021], [66.00222110470447, 5.617045168747429e-16, 9.173330910852394], [66.04640041397415, 5.510152142709861e-16, 8.998761351511707], [0.0, -0.16812879171659068, 1.0294919331011762e-17], 0.0)])
Alert!

Note that the fuselage is not included in the analysis for now, and is only included for plotting. Support for its aerodynamic analysis will be added soon.

Now, let's define the freestream conditions.

fs  = Freestream(
    alpha = 3.0, # degrees
    beta  = 0.0, # degrees
    omega = [0., 0., 0.]
);

Similarly, define the reference values. Here, the reference flight condition will be set to Mach number $M = 0.84$.

M = 0.84 # Mach number
refs = References(
    sound_speed = 330.,
    speed    = M * 330.,
    density  = 1.225,
    span     = b_w,
    area     = S_w,
    chord    = c_w,
    location = mac_w
);

Let's run the aerodynamic analysis first.

system = solve_case(
    aircraft, fs, refs;
    compressible     = true, # Compressibility option
    # print            = true, # Prints the results for only the aircraft
    # print_components = true, # Prints the results for all components
)
VortexLatticeSystem -
368 Horseshoe{Float64} Elements
Freestream: 
    alpha = 0.05235987755982989
    beta = 0.0
    omega = [0.0, 0.0, 0.0]
References: 
    Mach = 0.84
    Reynolds = 1.9693498716934252e8
    speed = 277.2
    density = 1.225
    viscosity = 1.5e-5
    sound_speed = 330.0
    area = 427.9435545
    span = 60.9
    chord = 8.699310326413222
    location = [30.683074339989055, 0.0, -2.5]
Info

You may receive a warning that the results are incorrect; this is due to the limitation of the vortex lattice method being able to primarily analyze only subsonic flows under the physical assumptions.

The freestream derivatives can be obtained by passing the resultant system as follows:

dvs = freestream_derivatives(
    system,                     # VortexLatticeSystem
    axes             = Wind(),  # Specify axis system for nearfield forces (wind by default)
    # print            = true,    # Prints the results for only the aircraft
    # print_components = true,    # Prints the results for all components
    # farfield         = true,    # Print farfield derivatives
);
┌ Warning: Results in transonic to sonic flow conditions (0.7 < M < 1) are most likely incorrect!
└ @ AeroFuse.VortexLattice ~/work/AeroFuse.jl/AeroFuse.jl/src/Aerodynamics/VortexLattice/system.jl:112

You can access the derivatives of each lifting surface based on the keys defined in the ComponentVector.

ac_dvs = dvs.aircraft
9×7 LabelledArrays.SLArray{Tuple{9, 7}, Float64, 2, 63, (:CX, :CY, :CZ, :Cl, :Cm, :Cn, :CDi, :CY, :CL, :CX_Ma, :CY_Ma, :CZ_Ma, :Cl_Ma, :Cm_Ma, :Cn_Ma, :CDi_Ma, :CY_Ma, :CL_Ma, :CX_al, :CY_al, :CZ_al, :Cl_al, :Cm_al, :Cn_al, :CDi_al, :CY_al, :CL_al, :CX_be, :CY_be, :CZ_be, :Cl_be, :Cm_be, :Cn_be, :CDi_be, :CY_be, :CL_be, :CX_pb, :CY_pb, :CZ_pb, :Cl_pb, :Cm_pb, :Cn_pb, :CDi_pb, :CY_pb, :CL_pb, :CX_qb, :CY_qb, :CZ_qb, :Cl_qb, :Cm_qb, :Cn_qb, :CDi_qb, :CY_qb, :CL_qb, :CX_rb, :CY_rb, :CZ_rb, :Cl_rb, :Cm_rb, :Cn_rb, :CDi_rb, :CY_rb, :CL_rb)} with indices SOneTo(9)×SOneTo(7):
  :CX => 0.0022641332528671202    …   :CX_rb => -1.5827149500136481e-12
  :CY => -3.1961713265252737e-16      :CY_rb => -3.157862528537969
  :CZ => 0.488779300272069            :CZ_rb => 7.974628825637478e-13
  :Cl => -3.250612744503939e-17       :Cl_rb => 0.33851067064045093
  :Cm => 0.01236035958459586          :Cm_rb => -3.074932469430325e-12
  :Cn => 1.6465191416848013e-16   …   :Cn_rb => 1.7043644190221912
 :CDi => 0.009363051346220326        :CDi_rb => 1.7542638011755205e-13
  :CY => -3.196582110622816e-16       :CY_rb => -3.153938658278606
  :CL => 0.4887854408983533           :CL_rb => 6.327189013599491e-15

These quantities are the force and moment coefficients $(C_X, C_Y, C_Z, C_l, C_m, C_n, C_{D_{i,ff}}, C_{Y_{ff}} C_{L_{ff}})$ generated from the nearfield and farfield analyses, and their derivatives respect to the Mach number $M$, freestream angles of attack and sideslip $(\alpha, \beta)$, and the non-dimensional angular velocity rates in stability axes $(\bar{p}, \bar{q}, \bar{r})$. The keys corresponding to the freestream derivatives should be evident:

keys(dvs.aircraft)
(:CX, :CY, :CZ, :Cl, :Cm, :Cn, :CDi, :CY, :CL, :CX_Ma, :CY_Ma, :CZ_Ma, :Cl_Ma, :Cm_Ma, :Cn_Ma, :CDi_Ma, :CY_Ma, :CL_Ma, :CX_al, :CY_al, :CZ_al, :Cl_al, :Cm_al, :Cn_al, :CDi_al, :CY_al, :CL_al, :CX_be, :CY_be, :CZ_be, :Cl_be, :Cm_be, :Cn_be, :CDi_be, :CY_be, :CL_be, :CX_pb, :CY_pb, :CZ_pb, :Cl_pb, :Cm_pb, :Cn_pb, :CDi_pb, :CY_pb, :CL_pb, :CX_qb, :CY_qb, :CZ_qb, :Cl_qb, :Cm_qb, :Cn_qb, :CDi_qb, :CY_qb, :CL_qb, :CX_rb, :CY_rb, :CZ_rb, :Cl_rb, :Cm_rb, :Cn_rb, :CDi_rb, :CY_rb, :CL_rb)

These can be accessed either like a dictionary, or by 'dot' syntax.

ac_dvs[:CZ_al], ac_dvs.CZ_al, ac_dvs.CL_al # Lift coefficient derivative wrt. alpha
(6.8101016101310154, 6.8101016101310154, 6.801125646406997)

Note that the nearfield forces and moments $(C_X, C_Y, C_Z, C_l, C_m, C_n)$ depend on the axis system used ($C_Z$ is not lift if body axes are used!). You can also pretty-print the derivatives for each surface.

print_derivatives(dvs.aircraft, "Aircraft", farfield = true)
print_derivatives(dvs.wing, "Wing", farfield = true)
print_derivatives(dvs.htail, "Horizontal Tail", farfield = true)
print_derivatives(dvs.vtail, "Vertical Tail", farfield = true)
─────────────────────────────────────────────────────────────────────────────────────────────
 Aircraft    Values                             Freestream   Derivatives
                          ∂/∂M     ∂/∂α, 1/rad  ∂/∂β, 1/rad     ∂/∂p̄        ∂/∂q̄      ∂/∂r̄
─────────────────────────────────────────────────────────────────────────────────────────────
    CX     0.00226413  0.00627624   0.0870775       0.0         -0.0      -1.69412    -0.0
    CY        -0.0        -0.0         0.0        2.61907     -0.556502     0.0     -3.15786
    CZ      0.488779    0.637111     6.8101        -0.0          0.0      113.021     0.0
    Cℓ        -0.0        0.0          0.0      -0.0662498    -0.524179     0.0     0.338511
    Cm     0.0123604   -0.0528664   -2.82087        0.0         -0.0      -92.7847    -0.0
    Cn        0.0         0.0         -0.0       -1.38315     0.224446      -0.0    1.70436
   CDi     0.00936305  0.0246163    0.195791       -0.0          0.0      2.55536     0.0
    CY        -0.0        -0.0         0.0        2.60569     -0.658042     0.0     -3.15394
    CL      0.488785    0.636387     6.80113        0.0          0.0      112.523     0.0
─────────────────────────────────────────────────────────────────────────────────────────────
────────────────────────────────────────────────────────────────────────────────────────────
 Wing    Values                             Freestream   Derivatives
                      ∂/∂M     ∂/∂α, 1/rad  ∂/∂β, 1/rad     ∂/∂p̄        ∂/∂q̄       ∂/∂r̄
────────────────────────────────────────────────────────────────────────────────────────────
  CX   0.00334069  0.00925357   0.0763912      -0.0         -0.0      -2.68788      0.0
  CY      -0.0        -0.0         0.0      -0.0402579   -0.0387053     0.0      0.113127
  CZ    0.52069     0.706475     6.06875        0.0          0.0      88.5232      -0.0
  Cℓ      -0.0        0.0         -0.0       -0.132842    -0.50224      -0.0     0.387917
  Cm   -0.0960197  -0.290382    -0.311637      -0.0          0.0      -9.07936     -0.0
  Cn      -0.0        -0.0        -0.0      -0.0021049   -0.0452094     -0.0    -0.00205007
 CDi   0.00909333   0.023449    0.208213        0.0          0.0      2.96871      -0.0
  CY      -0.0        0.0         -0.0      -0.0341629    -0.145646     0.0      0.0852901
  CL    0.52054     0.705389     6.06549        0.0          0.0      88.1339      -0.0
────────────────────────────────────────────────────────────────────────────────────────────
────────────────────────────────────────────────────────────────────────────────────────────────────────
 Horizontal Tail    Values                               Freestream   Derivatives
                                  ∂/∂M      ∂/∂α, 1/rad  ∂/∂β, 1/rad     ∂/∂p̄        ∂/∂q̄       ∂/∂r̄
────────────────────────────────────────────────────────────────────────────────────────────────────────
       CX         -0.00107656  -0.00297734   0.0106863       0.0         -0.0      0.993763     -0.0
       CY             0.0          0.0         -0.0      -0.0862034    0.0109909     -0.0     0.130153
       CZ         -0.0319106   -0.0693637    0.741355        0.0         -0.0       24.4973      0.0
       Cℓ             0.0          0.0         -0.0      -0.0547872   0.00257755     -0.0     0.0569628
       Cm           0.10838     0.237515     -2.50923       -0.0         -0.0      -83.7053     -0.0
       Cn            -0.0         -0.0          0.0       0.0454955   -0.00568786     0.0     -0.068279
       CDi        0.00026972   0.00116723    -0.012422      -0.0          0.0      -0.413348     0.0
       CY             0.0          0.0         -0.0      -0.0691799   0.00944901     -0.0      0.11129
       CL         -0.0317545   -0.0690021    0.735637        0.0         -0.0       24.3893      0.0
────────────────────────────────────────────────────────────────────────────────────────────────────────
─────────────────────────────────────────────────────────────────────────────────────
 Vertical Tail  Values                     Freestream   Derivatives
                        ∂/∂M  ∂/∂α, 1/rad  ∂/∂β, 1/rad     ∂/∂p̄      ∂/∂q̄    ∂/∂r̄
─────────────────────────────────────────────────────────────────────────────────────
      CX         -0.0   -0.0      0.0          0.0         -0.0      0.0     -0.0
      CY         -0.0   -0.0      0.0        2.74554     -0.528787   0.0   -3.40114
      CZ         0.0    0.0      -0.0         -0.0          0.0      -0.0     0.0
      Cℓ         -0.0   -0.0      0.0        0.12138    -0.0245174   0.0   -0.106369
      Cm         -0.0   -0.0      0.0          0.0         -0.0      0.0     -0.0
      Cn         0.0    0.0      -0.0       -1.42654     0.275343    -0.0   1.77469
      CDi        0.0    0.0      -0.0         -0.0          0.0      -0.0     0.0
      CY         -0.0   -0.0      0.0        2.70903     -0.521845   0.0   -3.35052
      CL         0.0    0.0      -0.0         -0.0          0.0      -0.0     0.0
─────────────────────────────────────────────────────────────────────────────────────

Static Stability Analysis

You can evaluate the static stability of the aircraft using the various quantities computed in this process, following standard computations and definitions in aircraft design and stability analysis.

l_h = x_h - x_w                 # Horizontal tail moment arm
V_h = S_h / S_w * l_h / c_w     # Horizontal tail volume coefficient

l_v = x_v - x_w                 # Vertical tail moment arm
V_v = S_v / S_w * l_v / b_w     # Vertical tail volume coefficient

x_cp = -refs.chord * ac_dvs.Cm / ac_dvs.CZ # Center of pressure

x_np_lon = -refs.chord * ac_dvs.Cm_al / ac_dvs.CZ_al # Neutral point

# Position vectors
x_np, y_np, z_np = r_np = refs.location + [ x_np_lon; zeros(2) ]  # Neutral point
x_cp, y_cp, z_cp = r_cp = refs.location + [ x_cp; zeros(2) ]  # Center of pressure

@info "Horizontal TVC                V_h:" V_h
@info "Vertical TVC                  V_v:" V_v
@info "Wing Aerodynamic Center  x_ac (m):" x_w
@info "Neutral Point            x_np (m):" x_np
@info "Center of Pressure       x_cp (m):" x_cp
┌ Info: Horizontal TVC                V_h:
└   V_h = 0.7977744718731725
┌ Info: Vertical TVC                  V_v:
└   V_v = 0.06132559058467538
┌ Info: Wing Aerodynamic Center  x_ac (m):
└   x_w = 30.683074339989055
┌ Info: Neutral Point            x_np (m):
└   x_np = 34.286491892544305
┌ Info: Center of Pressure       x_cp (m):
└   x_cp = 30.463084246888542

Let's plot the relevant quantities.

stab_plt = plot(
    xaxis = L"x", yaxis = L"y", zaxis = L"z",
    aspect_ratio = 1,
    zlim = (-0.5, 0.5) .* span(wing_mesh),
    camera = (30,60),
)
plot!(fuse, label = "Fuselage", alpha = 0.6)
plot!(stab_plt, wing_mesh, label = "Wing", mac = false)
plot!(stab_plt, htail_mesh, label = "Horizontal Tail", mac = false)
plot!(stab_plt, vtail_mesh, label = "Vertical Tail", mac = false)

scatter!(Tuple(r_np), color = :orange, label = "Neutral Point")
scatter!(Tuple(r_cp), color = :brown, label = "Center of Pressure")

This page was generated using Literate.jl.