Objectives
Here we will show you how to perform an aerodynamic analysis of a conventional aircraft.
Recipe
- Define the geometries of a wing, horizontal tail and vertical tail.
- Mesh and plot the geometries for numerical analyses.
- Perform an aerodynamic analysis of this aircraft configuration at given freestream conditions and reference values.
- Evaluate its drag polar for a given range of angles of attack.
- Plot the spanwise loading distribution for the wing.
For this, we will need to import some packages which will be convenient for plotting.
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
Your First Wing
Here you will learn how to define a wing using an intuitive parametrization scheme. First, we define a Vector
of Foil
s.
# Define one airfoil
airfoil_1 = naca4(4,4,1,2)
# Define vector of airfoils
airfoils = [ airfoil_1, naca4((0,0,1,2)) ]
2-element Vector{Foil{Float64}}:
NACA 4412 Foil{Float64} with 79 points.
NACA 0012 Foil{Float64} with 79 points.
Refer to the Airfoil Aerodynamic Analysis tutorial for an introduction to the Foil
type.
Parametrization
The following function defines a symmetric, single-section wing and prints the relevant information.
All relevant units are presented in metric.
wing = Wing(
foils = airfoils, # Foil profiles
chords = [1.0, 0.6], # Chord lengths
twists = [2.0, 0.0], # Twist angles (degrees)
spans = [4.0], # Section span lengths
dihedrals = [5.], # Dihedral angles (degrees)
sweeps = [5.], # Sweep angles (degrees)
w_sweep = 0., # Sweep angle location w.r.t.
# normalized chord lengths ∈ [0,1]
symmetry = true, # Whether wing is symmetric
# flip = false # Whether wing is reflected
)
Symmetric AbstractWing with 1 spanwise section(s).
Aspect Ratio: 10.001523280439075
Span (m): 8.0
Projected Area (m²): 6.399025249000904
Mean Aerodynamic Chord (m): 0.8166666666666667
Mean Aerodynamic Center (m): [0.36456254979752734, 0.0, 0.0]
Position (m): [0.0, 0.0, 0.0]
In this case, the root chord is $1.0~m$ and the tip chord is $0.6~m$, and so on.
See the How-to Guide for the various ways of constructing wings.
Visualization
You can plot the lifting surface by calling plot
with Plots.jl
.
plt = plot(
xlabel = L"x", ylabel = L"y", zlabel = L"z",
aspect_ratio = 1,
camera = (30, 30),
zlim = (-0.5, 0.5) .* span(wing),
)
plot!(plt, wing, label = "Wing")
Your First Vortex Lattice Analysis
Now we would like to analyze the aerodynamics of this wing in conjunction with other lifting surfaces pertinent to aircraft.
Geometry
We define the horizontal tail similarly to the wing. However, we also add additional position (by specifying a vector) and orientation attributes (by specifying an angle and axis of rotation) to place it at the desired location.
# Horizontal tail
htail = Wing(
foils = fill(naca4(0,0,1,2), 2),
chords = [0.7, 0.42],
twists = [0.0, 0.0],
spans = [1.25],
dihedrals = [0.],
sweeps = [6.39],
w_sweep = 0.,
position = [4., 0, 0],
angle = -2.,
axis = [0., 1., 0.],
symmetry = true
)
Symmetric AbstractWing with 1 spanwise section(s).
Aspect Ratio: 4.464285714285714
Span (m): 2.5
Projected Area (m²): 1.4
Mean Aerodynamic Chord (m): 0.5716666666666665
Mean Aerodynamic Center (m): [4.206952171902948, 0.0, 0.007226929090839386]
Position (m): [4.0, 0.0, 0.0]
For the vertical tail, we simply set symmetry = false
to define its shape.
# Vertical tail
vtail = Wing(
foils = fill(naca4(0,0,0,9), 2),
chords = [0.7, 0.42],
twists = [0.0, 0.0],
spans = [1.0],
dihedrals = [0.],
sweeps = [7.97],
w_sweep = 0.,
position = [4., 0, 0],
angle = 90.,
axis = [1., 0., 0.]
)
AbstractWing with 1 spanwise section(s).
Aspect Ratio: 1.7857142857142858
Span (m): 1.0
Projected Area (m²): 0.5599999999999999
Mean Aerodynamic Chord (m): 0.5716666666666665
Mean Aerodynamic Center (m): [4.207086511043837, 2.806482248046018e-17, 0.45833333333333337]
Position (m): [4.0, 0.0, 0.0]
Let's visualize the geometry of the aircraft's configuration.
plot!(plt, htail, label = "Horizontal Tail")
plot!(plt, vtail, label = "Vertical Tail")
Meshing
To perform the aerodynamic analysis, we will need to discretize our geometry into a mesh. The following WingMesh
function constructs a mesh for you by providing a HalfWing
or Wing
type with specification of spanwise panels and chordwise panels. As the wing has only one spanwise section, we provide a vector with a single integer entry for the spanwise panel distribution.
wing_mesh = WingMesh(wing, [20], 12) # (Wing, [Spanwise panels], Chordwise panels)
WingMesh —
Symmetric AbstractWing with 1 spanwise section(s).
Aspect Ratio: 10.001523280439075
Span (m): 8.0
Projected Area (m²): 6.399025249000904
Mean Aerodynamic Chord (m): 0.8166666666666667
Mean Aerodynamic Center (m): [0.36456254979752734, 0.0, 0.0]
Position (m): [0.0, 0.0, 0.0]
Spanwise panels: [10, 10]
Chordwise panels: 12
Spanwise spacing: AbstractSpacing[Sine(true), Sine(false)]
Chordwise spacing: Cosine()
Let's see what this discretization looks like on the camber distribution of the wing.
# Compute camber panel distribution
plot!(plt, wing_mesh, label = "")
Similarly we define the meshes for the other surfaces and plot them.
htail_mesh = WingMesh(htail, [8], 6)
vtail_mesh = WingMesh(vtail, [6], 4)
plot!(plt, htail_mesh, label = "")
plot!(plt, vtail_mesh, label = "")
Aerodynamic Analysis
For the analysis, we need to generate a Horseshoe
type, corresponding to horseshoe singularity elements used in the vortex lattice method. This is generated as follows:
make_horseshoes(wing_mesh)
12×20 Matrix{Horseshoe{Float64}}:
Horseshoe{Float64}([0.35251, -4.0, 0.349955], [0.348223, -3.95075, 0.345644], [0.355499, -3.97538, 0.347798], [-3.13285e-6, 8.75328e-5, 0.00101107], 0.0) … Horseshoe{Float64}([0.348223, 3.95075, 0.345644], [0.35251, 4.0, 0.349955], [0.355499, 3.97538, 0.347798], [-3.13285e-6, -8.75328e-5, 0.00101107], 0.0)
Horseshoe{Float64}([0.367669, -4.0, 0.349955], [0.363506, -3.95075, 0.345633], [0.380634, -3.97538, 0.347789], [-3.84966e-6, 0.000252825, 0.00296412], 0.0) Horseshoe{Float64}([0.363506, 3.95075, 0.345633], [0.367669, 4.0, 0.349955], [0.380634, 3.97538, 0.347789], [-3.84966e-6, -0.000252825, 0.00296412], 0.0)
Horseshoe{Float64}([0.402066, -4.0, 0.349955], [0.398185, -3.95075, 0.345609], [0.424061, -3.97538, 0.347773], [-4.85621e-6, 0.000393857, 0.00471519], 0.0) Horseshoe{Float64}([0.398185, 3.95075, 0.345609], [0.402066, 4.0, 0.349955], [0.424061, 3.97538, 0.347773], [-4.85621e-6, -0.000393857, 0.00471519], 0.0)
Horseshoe{Float64}([0.453356, -4.0, 0.349955], [0.449896, -3.95075, 0.345572], [0.482819, -3.97538, 0.347752], [-3.74655e-6, 0.000502611, 0.00614487], 0.0) Horseshoe{Float64}([0.449896, 3.95075, 0.345572], [0.453356, 4.0, 0.349955], [0.482819, 3.97538, 0.347752], [-3.74655e-6, -0.000502611, 0.00614487], 0.0)
Horseshoe{Float64}([0.518043, -4.0, 0.349955], [0.515114, -3.95075, 0.345526], [0.552904, -3.97538, 0.347727], [-6.19919e-7, 0.00057911, 0.00715572], 0.0) Horseshoe{Float64}([0.515114, 3.95075, 0.345526], [0.518043, 4.0, 0.349955], [0.552904, 3.97538, 0.347727], [-6.19919e-7, -0.00057911, 0.00715572], 0.0)
Horseshoe{Float64}([0.59172, -4.0, 0.349955], [0.589396, -3.95075, 0.345473], [0.62954, -3.97538, 0.3477], [3.20304e-6, 0.000625718, 0.0076789], 0.0) … Horseshoe{Float64}([0.589396, 3.95075, 0.345473], [0.59172, 4.0, 0.349955], [0.62954, 3.97538, 0.3477], [3.20304e-6, -0.000625718, 0.0076789], 0.0)
Horseshoe{Float64}([0.669366, -4.0, 0.349955], [0.667679, -3.95075, 0.345417], [0.707505, -3.97538, 0.347672], [5.53088e-6, 0.000639568, 0.00767882], 0.0) Horseshoe{Float64}([0.667679, 3.95075, 0.345417], [0.669366, 4.0, 0.349955], [0.707505, 3.97538, 0.347672], [5.53088e-6, -0.000639568, 0.00767882], 0.0)
Horseshoe{Float64}([0.745689, -4.0, 0.349955], [0.744629, -3.95075, 0.345363], [0.781484, -3.97538, 0.347646], [7.20081e-6, 0.000614732, 0.00715545], 0.0) Horseshoe{Float64}([0.744629, 3.95075, 0.345363], [0.745689, 4.0, 0.349955], [0.781484, 3.97538, 0.347646], [7.20081e-6, -0.000614732, 0.00715545], 0.0)
Horseshoe{Float64}([0.815488, -4.0, 0.349955], [0.815, -3.95075, 0.345313], [0.846437, -3.97538, 0.347623], [7.80914e-6, 0.000546835, 0.00614447], 0.0) Horseshoe{Float64}([0.815, 3.95075, 0.345313], [0.815488, 4.0, 0.349955], [0.846437, 3.97538, 0.347623], [7.80914e-6, -0.000546835, 0.00614447], 0.0)
Horseshoe{Float64}([0.874006, -4.0, 0.349955], [0.873998, -3.95075, 0.345271], [0.897938, -3.97538, 0.347604], [7.01154e-6, 0.000433948, 0.00471478], 0.0) Horseshoe{Float64}([0.873998, 3.95075, 0.345271], [0.874006, 4.0, 0.349955], [0.897938, 3.97538, 0.347604], [7.01154e-6, -0.000433948, 0.00471478], 0.0)
Horseshoe{Float64}([0.917255, -4.0, 0.349955], [0.917603, -3.95075, 0.34524], [0.932475, -3.97538, 0.347592], [4.86651e-6, 0.000280016, 0.00296382], 0.0) … Horseshoe{Float64}([0.917603, 3.95075, 0.34524], [0.917255, 4.0, 0.349955], [0.932475, 3.97538, 0.347592], [4.86651e-6, -0.000280016, 0.00296382], 0.0)
Horseshoe{Float64}([0.942288, -4.0, 0.349955], [0.942841, -3.95075, 0.345222], [0.947697, -3.97538, 0.347586], [1.7431e-6, 9.68789e-5, 0.0010109], 0.0) Horseshoe{Float64}([0.942841, 3.95075, 0.345222], [0.942288, 4.0, 0.349955], [0.947697, 3.97538, 0.347586], [1.7431e-6, -9.68789e-5, 0.0010109], 0.0)
To perform the aerodynamic analysis, you have to assemble these horseshoes for each surface into a ComponentVector
.
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}([0.3525102171320159, -4.0, 0.34995465410369603], [0.3482226701046064, -3.950753362380551, 0.34564430151089465], [0.3554985451560125, -3.9753766811902755, 0.34779764772296495], [-3.132846635733752e-6, 8.75328405578299e-5, 0.0010110706764502092], 0.0) Horseshoe{Float64}([0.3482226701046064, -3.950753362380551, 0.34564430151089465], [0.33546560265974534, -3.804226065180614, 0.3328193789185651], [0.34705962335816604, -3.877489713780583, 0.339222734855828], [-4.63774703535055e-5, 0.00026159674585068623, 0.00305855758747353], 0.0) … Horseshoe{Float64}([0.33546560265974534, 3.804226065180614, 0.3328193789185651], [0.3482226701046064, 3.950753362380551, 0.34564430151089465], [0.34705962335816604, 3.877489713780583, 0.339222734855828], [-4.63774703535055e-5, -0.00026159674585068623, 0.00305855758747353], 0.0) Horseshoe{Float64}([0.3482226701046064, 3.950753362380551, 0.34564430151089465], [0.3525102171320159, 4.0, 0.34995465410369603], [0.3554985451560125, 3.9753766811902755, 0.34779764772296495], [-3.132846635733752e-6, -8.75328405578299e-5, 0.0010110706764502092], 0.0); Horseshoe{Float64}([0.36766943790482276, -4.0, 0.34995465410369603], [0.3635063143194118, -3.950753362380551, 0.34563344572213517], [0.38063443648662565, -3.9753766811902755, 0.3477886843771474], [-3.8496592168059765e-6, 0.00025282532658469487, 0.0029641229179659502], 0.0) Horseshoe{Float64}([0.3635063143194118, -3.950753362380551, 0.34563344572213517], [0.35111945348246687, -3.804226065180614, 0.3327762230690753], [0.37260391801049586, -3.877489713780583, 0.33917813883406156], [-5.698889124096818e-5, 0.0007608534505138116, 0.008964478971185628], 0.0) … Horseshoe{Float64}([0.35111945348246687, 3.804226065180614, 0.3327762230690753], [0.3635063143194118, 3.950753362380551, 0.34563344572213517], [0.37260391801049586, 3.877489713780583, 0.33917813883406156], [-5.698889124096818e-5, -0.0007608534505138116, 0.008964478971185628], 0.0) Horseshoe{Float64}([0.3635063143194118, 3.950753362380551, 0.34563344572213517], [0.36766943790482276, 4.0, 0.34995465410369603], [0.38063443648662565, 3.9753766811902755, 0.3477886843771474], [-3.8496592168059765e-6, -0.00025282532658469487, 0.0029641229179659502], 0.0); … ; Horseshoe{Float64}([0.9172548069268747, -4.0, 0.34995465410369603], [0.9176025552596606, -3.950753362380551, 0.345239877823049], [0.9324752414677762, -3.9753766811902755, 0.34759190042760424], [4.866511340696338e-6, 0.0002800156739994455, 0.0029638185425830753], 0.0) Horseshoe{Float64}([0.9176025552596606, -3.950753362380551, 0.345239877823049], [0.9186372375399541, -3.804226065180614, 0.3312116424205826], [0.9334109305093639, -3.877489713780583, 0.3381990645602721], [7.204198343249224e-5, 0.0008460496495163087, 0.008959973113753865], 0.0) … Horseshoe{Float64}([0.9186372375399541, 3.804226065180614, 0.3312116424205826], [0.9176025552596606, 3.950753362380551, 0.345239877823049], [0.9334109305093639, 3.877489713780583, 0.3381990645602721], [7.204198343249224e-5, -0.0008460496495163087, 0.008959973113753865], 0.0) Horseshoe{Float64}([0.9176025552596606, 3.950753362380551, 0.345239877823049], [0.9172548069268747, 4.0, 0.34995465410369603], [0.9324752414677762, 3.9753766811902755, 0.34759190042760424], [4.866511340696338e-6, -0.0002800156739994455, 0.0029638185425830753], 0.0); Horseshoe{Float64}([0.9422879650187365, -4.0, 0.34995465410369603], [0.9428411798290255, -3.950753362380551, 0.34522195113141374], [0.9476966739615824, -3.9753766811902755, 0.3475864725332245], [1.7431038793475474e-6, 9.687894049861681e-5, 0.0010109004045062244], 0.0) Horseshoe{Float64}([0.9428411798290255, -3.950753362380551, 0.34522195113141374], [0.9444872022752491, -3.804226065180614, 0.3311403770686849], [0.9488796780281274, -3.877489713780583, 0.3381720587411475], [2.5804247027441836e-5, 0.00029263665761506297, 0.0030560369463593856], 0.0) … Horseshoe{Float64}([0.9444872022752491, 3.804226065180614, 0.3311403770686849], [0.9428411798290255, 3.950753362380551, 0.34522195113141374], [0.9488796780281274, 3.877489713780583, 0.3381720587411475], [2.5804247027441836e-5, -0.00029263665761506297, 0.0030560369463593856], 0.0) Horseshoe{Float64}([0.9428411798290255, 3.950753362380551, 0.34522195113141374], [0.9422879650187365, 4.0, 0.34995465410369603], [0.9476966739615824, 3.9753766811902755, 0.3475864725332245], [1.7431038793475474e-6, -9.687894049861681e-5, 0.0010109004045062244], 0.0)], htail = Horseshoe{Float64}[Horseshoe{Float64}([4.146933162399598, -1.25, 0.005131019094849938], [4.136640341075092, -1.1548494156391085, 0.004771585853957062], [4.156202234775822, -1.2024247078195542, 0.005454702234702354], [-0.0001915951254382746, 1.0164395367051604e-19, 0.00548656653982245], 0.0) Horseshoe{Float64}([4.136640341075092, -1.1548494156391085, 0.004771585853957062], [4.107328865843567, -0.8838834764831843, 0.0037480065839339315], [4.1377726586296575, -1.0193664460611465, 0.004811127254271502], [-0.0005975677270916761, 4.1504614415460717e-20, 0.017112100786693396], 0.0) … Horseshoe{Float64}([4.107328865843567, 0.8838834764831843, 0.0037480065839339315], [4.136640341075092, 1.1548494156391085, 0.004771585853957062], [4.1377726586296575, 1.0193664460611465, 0.004811127254271502], [-0.0005975677270916761, -4.1504614415460717e-20, 0.017112100786693396], 0.0) Horseshoe{Float64}([4.136640341075092, 1.1548494156391085, 0.004771585853957062], [4.146933162399598, 1.25, 0.005131019094849938], [4.156202234775822, 1.2024247078195542, 0.005454702234702354], [-0.0001915951254382746, -1.0164395367051604e-19, 0.00548656653982245], 0.0); Horseshoe{Float64}([4.18722593477255, -1.2500000000000002, 0.006538073711069197], [4.178977849894339, -1.1548494156391087, 0.0062500442402888285], [4.222485724410211, -1.2024247078195545, 0.007769372697333474], [-0.0005234476171798984, -4.0657581468206416e-19, 0.014989578545902158], 0.0) Horseshoe{Float64}([4.178977849894339, -1.1548494156391087, 0.0062500442402888285], [4.155489291413338, -0.8838834764831844, 0.005429805703880379], [4.210367339531874, -1.0193664460611465, 0.007346189372384801], [-0.0016325853913779116, -5.963111948670274e-19, 0.04675112877348567], 0.0) … Horseshoe{Float64}([4.155489291413338, 0.8838834764831844, 0.005429805703880379], [4.178977849894339, 1.1548494156391087, 0.0062500442402888285], [4.210367339531874, 1.0193664460611465, 0.007346189372384801], [-0.0016325853913779116, 5.963111948670274e-19, 0.04675112877348567], 0.0) Horseshoe{Float64}([4.178977849894339, 1.1548494156391087, 0.0062500442402888285], [4.18722593477255, 1.2500000000000002, 0.006538073711069197], [4.222485724410211, 1.2024247078195545, 0.007769372697333474], [-0.0005234476171798984, 4.0657581468206416e-19, 0.014989578545902158], 0.0); … ; Horseshoe{Float64}([4.473916518956167, -1.25, 0.016549529516799787], [4.480217114939736, -1.1548494156391085, 0.01676955117680262], [4.516450649024717, -1.2024247078195542, 0.018034854068455668], [-0.0005234476171798977, 2.4123498337802474e-18, 0.014989578545902158], 0.0) Horseshoe{Float64}([4.480217114939736, -1.1548494156391085, 0.01676955117680262], [4.498159694266988, -0.8838834764831843, 0.01739611985357701], [4.532322173481397, -1.0193664460611465, 0.01858909991549001], [-0.0016325853913779107, -3.415236843329339e-18, 0.0467511287734859], 0.0) … Horseshoe{Float64}([4.498159694266988, 0.8838834764831843, 0.01739611985357701], [4.480217114939736, 1.1548494156391085, 0.01676955117680262], [4.532322173481397, 1.0193664460611465, 0.01858909991549001], [-0.0016325853913779107, 3.415236843329339e-18, 0.0467511287734859], 0.0) Horseshoe{Float64}([4.480217114939736, 1.1548494156391085, 0.01676955117680262], [4.473916518956167, 1.25, 0.016549529516799787], [4.516450649024717, 1.2024247078195542, 0.018034854068455668], [-0.0005234476171798977, -2.4123498337802474e-18, 0.014989578545902158], 0.0); Horseshoe{Float64}([4.538559783420221, -1.25, 0.01880692205434313], [4.548140829744457, -1.1548494156391085, 0.019141499564521516], [4.557765789620817, -1.2024247078195542, 0.01947761056973118], [-0.00019159512543827543, 7.216720710606639e-19, 0.00548656653982245], 0.0) Horseshoe{Float64}([4.548140829744457, -1.1548494156391085, 0.019141499564521516], [4.5754253412667, -0.8838834764831843, 0.020094295702084708], [4.577571140675906, -1.0193664460611465, 0.02016922866862912], [-0.0005975677270916782, -9.825582188149884e-20, 0.017112100786693396], 0.0) … Horseshoe{Float64}([4.5754253412667, 0.8838834764831843, 0.020094295702084708], [4.548140829744457, 1.1548494156391085, 0.019141499564521516], [4.577571140675906, 1.0193664460611465, 0.02016922866862912], [-0.0005975677270916782, 9.825582188149884e-20, 0.017112100786693396], 0.0) Horseshoe{Float64}([4.548140829744457, 1.1548494156391085, 0.019141499564521516], [4.538559783420221, 1.25, 0.01880692205434313], [4.557765789620817, 1.2024247078195542, 0.01947761056973118], [-0.00019159512543827543, -7.216720710606639e-19, 0.00548656653982245], 0.0)], vtail = Horseshoe{Float64}[Horseshoe{Float64}([4.025628156646177, 0.0, 0.0], [4.034320138429017, 4.101789010561156e-18, 0.06698729810778065], [4.0805437564422675, 2.050894505280578e-18, 0.03349364905389032], [0.0, -0.013550085867564059, 8.2970346429420555e-19], 0.0) Horseshoe{Float64}([4.034320138429017, 4.101789010561156e-18, 0.06698729810778065], [4.058067074278198, 1.530808498934191e-17, 0.2499999999999999], [4.09420039959366, 9.704936999951532e-18, 0.15849364905389027], [0.0, -0.03514341176017239, 2.1519133361606285e-18], 0.0) … Horseshoe{Float64}([4.122944909542239, 4.592425496802574e-17, 0.75], [4.1466918453914206, 5.713055094680649e-17, 0.9330127018922192], [4.168821685494544, 5.1527402957416116e-17, 0.8415063509461096], [0.0, -0.024892149101701472, 1.5242045360648685e-18], 0.0) Horseshoe{Float64}([4.1466918453914206, 5.713055094680649e-17, 0.9330127018922192], [4.15538382717426, 6.123233995736766e-17, 1.0], [4.182478328645937, 5.918144545208708e-17, 0.9665063509461096], [0.0, -0.008424454538328645, 5.158490642463289e-19], 0.0); Horseshoe{Float64}([4.164384469938531, 0.0, 0.0], [4.1693584875122305, 4.101789010561156e-18, 0.06698729810778065], [4.288957314386939, 2.050894505280578e-18, 0.03349364905389032], [0.0, -0.03271280107279308, 2.0030813562470077e-18], 0.0) Horseshoe{Float64}([4.1693584875122305, 4.101789010561156e-18, 0.06698729810778065], [4.182947756241317, 1.530808498934191e-17, 0.2499999999999999], [4.292051773202949, 9.704936999951532e-18, 0.15849364905389027], [0.0, -0.0848437012994703, 5.195178361210522e-18], 0.0) … Horseshoe{Float64}([4.220074328846888, 4.592425496802574e-17, 0.75], [4.2336635975759735, 5.713055094680649e-17, 0.9330127018922192], [4.308960210617489, 5.1527402957416116e-17, 0.8415063509461096], [0.0, -0.060094963957941214, 3.679755262798413e-18], 0.0) Horseshoe{Float64}([4.2336635975759735, 5.713055094680649e-17, 0.9330127018922192], [4.238637615149672, 6.123233995736766e-17, 1.0], [4.312054669433499, 5.918144545208708e-17, 0.9665063509461096], [0.0, -0.020338432402028615, 1.245369807040959e-18], 0.0); Horseshoe{Float64}([4.411871843353823, 0.0, 0.0], [4.410214456743169, 4.101789010561156e-18, 0.06698729810778065], [4.533128985710054, 2.050894505280578e-18, 0.03349364905389032], [0.0, -0.03271280107279314, 2.0030813562470108e-18], 0.0) Horseshoe{Float64}([4.410214456743169, 4.101789010561156e-18, 0.06698729810778065], [4.405686392315079, 1.530808498934191e-17, 0.2499999999999999], [4.523849075855299, 9.704936999951532e-18, 0.15849364905389027], [0.0, -0.0848437012994703, 5.195178361210523e-18], 0.0) … Horseshoe{Float64}([4.393315490237592, 4.592425496802574e-17, 0.75], [4.388787425809501, 5.713055094680649e-17, 0.9330127018922192], [4.473142705429604, 5.1527402957416116e-17, 0.8415063509461096], [0.0, -0.060094963957941214, 3.679755262798413e-18], 0.0) Horseshoe{Float64}([4.388787425809501, 5.713055094680649e-17, 0.9330127018922192], [4.387130039198848, 6.123233995736766e-17, 1.0], [4.463862795574849, 5.918144545208708e-17, 0.9665063509461096], [0.0, -0.02033843240202856, 1.2453698070409554e-18], 0.0); Horseshoe{Float64}([4.623115530061469, 0.0, 0.0], [4.6157978859248665, 4.101789010561156e-18, 0.06698729810778065], [4.670026316897838, 2.050894505280578e-18, 0.03349364905389032], [0.0, -0.01355008586756406, 8.297034642942056e-19], 0.0) Horseshoe{Float64}([4.6157978859248665, 4.101789010561156e-18, 0.06698729810778065], [4.59580571035196, 1.530808498934191e-17, 0.2499999999999999], [4.653808591378466, 9.704936999951532e-18, 0.15849364905389027], [0.0, -0.035143411760172556, 2.1519133361606385e-18], 0.0) … Horseshoe{Float64}([4.541186070932944, 4.592425496802574e-17, 0.75], [4.521193895360037, 5.713055094680649e-17, 0.9330127018922192], [4.5651932911742055, 5.1527402957416116e-17, 0.8415063509461096], [0.0, -0.024892149101701472, 1.5242045360648685e-18], 0.0) Horseshoe{Float64}([4.521193895360037, 5.713055094680649e-17, 0.9330127018922192], [4.513876251223435, 6.123233995736766e-17, 1.0], [4.548975565654833, 5.918144545208708e-17, 0.9665063509461096], [0.0, -0.008424454538328586, 5.158490642463252e-19], 0.0)])
You can define the freestream condition as follows, by providing the angles of attack $\alpha$ and sideslip $\beta$ in degrees with a rotation vector $\Omega$.
fs = Freestream(
alpha = 3.0, # degrees
beta = 0.0, # degrees
omega = [0., 0., 0.]
);
You can define the reference values for the speed, area, span, chord, density, and location as follows.
refs = References(
speed = 1.0,
area = projected_area(wing),
span = span(wing),
chord = mean_aerodynamic_chord(wing),
density = 1.225,
location = mean_aerodynamic_center(wing)
);
You can run the aerodynamic analysis by providing the aircraft configuration, freestream, and reference values. Optionally you can also print the results.
system = solve_case(
aircraft, fs, refs;
compressible = false, # Compressibility option
# print = true, # Prints the results for only the aircraft
# print_components = true, # Prints the results for all components
)
VortexLatticeSystem -
312 Horseshoe{Float64} Elements
Freestream:
alpha = 0.05235987755982989
beta = 0.0
omega = [0.0, 0.0, 0.0]
References:
Mach = 0.0030303030303030303
Reynolds = 66694.44444444445
speed = 1.0
density = 1.225
viscosity = 1.5e-5
sound_speed = 330.0
area = 6.399025249000904
span = 8.0
chord = 0.8166666666666667
location = [0.36456254979752734, 0.0, 0.0]
You can obtain the aerodynamic coefficients from this system. The nearfield aerodynamic force and moment coefficients are ordered as $(C_{D_i}, C_Y, C_L, C_\ell, C_m, C_n)$.
nf = nearfield(system)
6-element LabelledArrays.SLArray{Tuple{6}, Float64, 1, 6, (:CX, :CY, :CZ, :Cl, :Cm, :Cn)} with indices SOneTo(6):
:CX => 0.009353793030993264
:CY => 1.1892283291375868e-15
:CZ => 0.5761421108553897
:Cl => -3.374056696043956e-17
:Cm => 0.07877036870774405
:Cn => -6.100017955052987e-16
The force coefficients are printed as $(C_X, C_Y, C_Z)$ for general axis systems; wind axes are used in the nearfield
function.
Refer to the how-to guide to see how to compute the aerodynamic coefficients of each component and perform stability analyses.
A convenience method is also provided for plotting streamlines from the leading edge of each surface.
plot!(plt,
system, # VortexLattice System
wing_mesh, # Lifting surface (or mesh)
span = 4, # Number of streamlines per spanwise section
dist = 10, # Distance of streamlines (m)
# lc = :green, # Color of streamlines
)
plot!(plt, system, htail_mesh, span = 3, lc = :cyan) # For horizontal tail
plot!(plt, system, vtail, span = 2, lc = :cyan) # For vertical tail
Drag Polar
Now let's analyze the drag polar of this aircraft configuration by varying the angle of attack and collecting the induced drag coefficient $C_{D_i}$.
# Define function to compute system varying with angle of attack.
vary_alpha(aircraft, α, refs) = VortexLatticeSystem(aircraft, Freestream(alpha = α), refs)
# Run loop
αs = -5:0.5:5
systems = [ vary_alpha(aircraft, α, refs) for α in αs ]
# Cleaner: map(α -> vary_alpha(...), αs)
# Get coefficients
coeffs = nearfield.(systems)
CDis = [ c[1] for c in coeffs ]
CLs = [ c[3] for c in coeffs ];
Let's plot the drag polar!
plot(CDis, CLs,
label = "",
xlabel = L"C_{D_i}",
ylabel = L"C_L",
title = "Drag Polar",
ls = :solid)
Let's also take a look at the variations of all the coefficients.
# Concatenate results into one array
data = permutedims(
mapreduce(hcat, systems) do sys
[sys.freestream.alpha; nearfield(sys) ]
end
)
# Plot
plot(
data[:,1], # Angle of attack
round.(data[:,2:end], digits = 4), # Aerodynamic coefficients
layout = (3,2),
xlabel = L"\alpha",
ylabel = [L"C_{D_i}" L"C_Y" L"C_L" L"C_\ell" L"C_m" L"C_n"],
labels = "",
)
Tip: You can convert this into a DataFrame for convenient reference.
using DataFrames
df = DataFrame(round.(data, digits = 6), [:α, :CX, :CY, :CZ, :Cl, :Cm, :Cn])
Spanwise Loading
You can compute the aerodynamic coefficients on the panels from the system.
CFs, CMs = surface_coefficients(system)
# Compute spanwise loads
span_loads = spanwise_loading(wing_mesh, system.reference, CFs.wing, system.circulations.wing)
# Plot spanwise loadings
plot_CD = plot(span_loads[:,1], span_loads[:,2], label = :none, ylabel = L"C_{D_i}")
plot_CY = plot(span_loads[:,1], span_loads[:,3], label = :none, ylabel = L"C_Y")
plot_CL = begin
plot(span_loads[:,1], span_loads[:,4], label = :none, xlabel = L"y", ylabel = L"C_L")
plot!(span_loads[:,1], span_loads[:,5], label = L"Γ/ρVc", xlabel = L"y")
end
plot(plot_CD, plot_CY, plot_CL, size = (800, 700), layout = (3,1))
This page was generated using Literate.jl.