How-to Guide

Airfoil Geometry

How to work with airfoil geometry.

Import Coordinates File

You can specify the path consisting of the foil's coordinates to the read_foil function. The Selig format for the coordinates file is followed, in which a header for the name in the first line of the file with the coordinates following from the second line are required. A custom name can be provided by setting the optional name variable.

# Airfoil coordinates file path
foilpath = string(@__DIR__, "/misc/s1223.dat")
"/home/runner/work/AeroFuse.jl/AeroFuse.jl/docs/build/misc/s1223.dat"

Note that, in this example, @__DIR__ points to the Julia REPL's current working directory. Use the correct local path for the coordinates on your computer.

# Read coordinates file
my_foil = read_foil(
    foilpath; # Path to the airfoil coordinates file
    name = "S1223" # To overwrite name in header
)
S1223 Foil{Float16} with 300 points.
plot(xlabel = L"x", ylabel = L"y", aspect_ratio = 1)
plot!(my_foil)
Tip

You can also import airfoil coordinates from the internet.

# Read coordinates from an internet address.
online_foil = read_foil(download("https://m-selig.ae.illinois.edu/ads/coord/s1210.dat"))

plot!(online_foil)

Interpolate and Process Coordinates

A basic cosine interpolation functionality is provided for foils.

# Cosine spacing with approx. 51 points on upper and lower surfaces each
cos_foil = cosine_interpolation(my_foil, 51)
S1223 Foil{Float64} with 99 points.

The upper and lower surfaces can be obtained by the following variety of functions.

# Upper surface
upper = upper_surface(my_foil)

# Lower surface
lower = lower_surface(my_foil)

# Upper and lower surfaces
upper, lower = split_surface(my_foil)
(Float16[-2.0e-5 -0.00073; -1.0e-5 0.00056; … ; 0.9976 0.00201; 1.0 0.0], Float16[-2.0e-5 -0.00073; 3.0e-5 -0.00202; … ; 0.997 0.00181; 1.0 0.0])

The camber-thickness distribution can be obtained as follows. It internally performs a cosine interpolation to standardize the $x$–coordinates for both surfaces; hence the number of points for the interpolation can be specified.

xcamthick = camber_thickness(my_foil, 60)
61×3 Matrix{Float64}:
 0.0          -0.00124659   0.0
 0.000685233  -9.97546e-5   0.0133951
 0.00273905    0.00142531   0.0254514
 0.00615583    0.00429351   0.0358441
 0.0109262     0.00779593   0.0455249
 0.0170371     0.0117405    0.0544609
 0.0244717     0.0157328    0.0629796
 0.0332098     0.0197507    0.0711215
 0.0432273     0.0239027    0.078745
 0.0544967     0.0281499    0.0859677
 ⋮                          
 0.956773      0.0258329    0.00703525
 0.96679       0.0215117    0.00653011
 0.975528      0.0171055    0.00556079
 0.982963      0.0126722    0.004475
 0.989074      0.00837342   0.00283016
 0.993844      0.00473417   0.0014028
 0.997261      0.00198852   0.000592446
 0.999315      0.000493805  0.000140883
 1.0           0.0          0.0

You can also do the inverse transformation.

coords = camber_thickness_to_coordinates(xcamthick[:,1], xcamthick[:,2], xcamthick[:,3])
121×2 Matrix{Float64}:
 1.0       0.0
 0.999315  0.000564246
 0.997261  0.00228474
 0.993844  0.00543557
 0.989074  0.0097885
 0.982963  0.0149097
 0.975528  0.0198859
 0.96679   0.0247768
 0.956773  0.0293506
 0.945503  0.0340995
 ⋮         
 0.956773  0.0223153
 0.96679   0.0182467
 0.975528  0.0143251
 0.982963  0.0104347
 0.989074  0.00695834
 0.993844  0.00403277
 0.997261  0.00169229
 0.999315  0.000423363
 1.0       0.0

Control Surfaces

You can (somewhat) mimic the behaviour of a control surface by specifying a deflection angle $\delta$ (in degrees, clockwise-positive convention) with the specification of the hinge location's $x$-coordinate normalized in $[0,1]$ to the chord length.

con_foil = control_surface(cos_foil; angle = 10., hinge = 0.75)

plot!(con_foil)

Doublet-Source Aerodynamic Analyses

The solve_case method runs the analysis given a Foil containing the airfoil coordinates, a Uniform2D defining the boundary conditions, and an optional named specification for the number of panels. It returns a system which can be used to obtain the aerodynamic quantities of interest and post-processing.

# Define freestream boundary conditions
uniform = Uniform2D(1.0, 4.0)

# Solve system
system  = solve_case(
    my_foil, # Foil
    uniform; # Freestream condition
    num_panels = 80
)
DoubletSourceSystem —
78 Panel2D{Float64} Elements
Freestream —
    alpha: 0.06981317007977318

The following functions compute the quantities of interest, such as the inviscid edge velocities, lift coefficient, and the sectional lift, moment, and pressure coefficients.

cls, cms, cps = surface_coefficients(system);
u_es   = surface_velocities(system)
cl     = lift_coefficient(system)
1.4921827407057007

AeroFuse provides more helper functions for the panel geometry.

panels   = system.surface_panels
pts      = collocation_point.(panels) # Collocation points
tangents = tangent_vector.(panels)     # Tangents
normals  = normal_vector.(panels)      # Normals
locs     = panel_location.(panels);   # Upper or lower surface

Wing Geometry

How to work with wing geometry.

To define a wing, AeroFuse provides a Wing constructor based on the following parametrization. The named arguments correspond to the foil shapes, chord and span lengths, twist, dihedral and sweep angles.

Info

A wing section consists of two foil profiles and their chord lengths and twist angles. Between them is their span length with associated leading-edge dihedral and sweep angles. So a general half-wing consisting of $n$ sections will have $n$ entries for spans $b$, dihedrals $\delta$, sweeps $\Lambda$, and $n+1$ entries for foils, chords $c$, and twists $\iota$, for some $n \in \mathbb N$.

airfoil = naca4((2,4,1,2))
wing = Wing(
    foils     = [ airfoil for i in 1:3 ],   # Foils
    chords    = [0.4, 0.2, 0.1],  # Chord lengths
    twists    = [0., 2., 5.],     # Twist angles (degrees)
    spans     = [1.0, 0.1],       # Section span lengths
    dihedrals = [0., 60.],        # Dihedral angles (degrees)
    sweeps    = [0., 30.],        # Sweep angles (degrees)
    w_sweep   = 0.25,             # Sweep angle location w.r.t.
                                  # normalized chord lengths ∈ [0,1]
    symmetry  = true,             # Whether wing is symmetric
    # flip      = false           # Whether wing is flipped in x-z plane
)
Symmetric AbstractWing with 2 spanwise section(s).
Aspect Ratio: 7.684336828804395
Span (m): 2.2
Projected Area (m²): 0.6298526610464908
Mean Aerodynamic Chord (m): 0.3037157904580297
Mean Aerodynamic Center (m): [0.10257127268682563, 0.0, 0.0]
Position (m): [0.0, 0.0, 0.0]

The symmetry Boolean argument specifies whether the geometry should be reflected in the $x$-$z$ plane. The flip Boolean argument specifies whether the coordinates should be flipped in the $x$-$z$ plane. The w_sweep argument specifies the chordwise-ratio of the sweep angles, e.g. 0. = leading edge sweep angle (default), 1. = trailing edge, 0.25 = quarter-chord.

The following "getter" functions provide quantities of interest such as chord lengths, spans, twist, dihedral, and sweep angles.

julia> foils(wing)3-element Vector{Foil{Float64}}:
 NACA 2412 Foil{Float64} with 79 points.
 NACA 2412 Foil{Float64} with 79 points.
 NACA 2412 Foil{Float64} with 79 points.
julia> chords(wing)3-element Vector{Float64}: 0.4 0.2 0.1
julia> spans(wing)2-element Vector{Float64}: 1.0 0.1
julia> twists(wing)3-element Vector{Float64}: -0.0 -2.0 -5.0
julia> dihedrals(wing)2-element Vector{Float64}: 0.0 60.0
julia> sweeps(wing) # Leading-edge sweep angles2-element Vector{Float64}: 2.862405226111748 44.03624346792648
julia> sweeps(wing, 0.25) # Quarter-chord sweep angles2-element Vector{Float64}: 0.0 30.000000000000004
julia> sweeps(wing, 1.0) # Trailing-edge sweep angles2-element Vector{Float64}: -8.447527247908468 -0.9637565320735177

You can evaluate commonly defined properties such as the aspect ratio, projected area, etc., with the following functions.

julia> aspect_ratio(wing)7.684336828804395
julia> projected_area(wing)0.6298526610464908
julia> span(wing)2.2
julia> taper_ratio(wing)0.25
julia> mean_aerodynamic_chord(wing)0.3037157904580297

The mean aerodynamic center calculation merits discussion. Its calculation procedure is depicted in the following diagram:

mac25_w = mean_aerodynamic_center(wing) # 25% MAC by default
mac40_w = mean_aerodynamic_center(wing, 0.4) # at 40% MAC
3-element StaticArraysCore.SVector{3, Float64} with indices SOneTo(3):
 0.14812864125553007
 0.0
 0.0

You can plot the wing with Plots.jl quite simply.

plot(wing,
    0.25; # Aerodynamic center location in terms of chord ratio plot
    mac = true, # Optional named argument to disable aerodynamic center plot
    zlim = (-0.5, 0.5) .* span(wing),
    aspect_ratio = 1,
    label = "Wing",
)

At 40% MAC

plot!(wing,
    0.40; # 40% MAC
    label = "40% MAC"
)

There is also a convenient function for pretty-printing information (using PrettyTables.jl), in which the first argument takes the Wing type and the second takes a name (as a String or Symbol).

print_info(wing, "Wing")
───────────────────────────────────────────────────
            Wing                     Value
───────────────────────────────────────────────────
        Aspect Ratio                7.68434
          Span (m)                    2.2
          Area (m²)                 0.629853
 Mean Aerodynamic Chord (m)         0.303716
 Mean Aerodynamic Center (m)  [0.102571, 0.0, 0.0]
───────────────────────────────────────────────────
Tip

You can use Accessors.jl to conveniently copy and modify properties of an existing object.

# Import Accessors
using Accessors

# Set only chords with other properties remaining identical.
new_wing = @set wing.chords = [0.4, 0.1, 0.05]
Symmetric AbstractWing with 2 spanwise section(s).
Aspect Ratio: 9.399958877486162
Span (m): 2.2
Projected Area (m²): 0.5148958695545236
Mean Aerodynamic Chord (m): 0.27411982939806356
Mean Aerodynamic Center (m): [0.09129804700674422, 0.0, 0.0]
Position (m): [0.0, 0.0, 0.0]

Vortex Lattice Aerodynamic Analyses

How to run a generic aerodynamic analysis on a conventional aircraft configuration.

Geometry

First we define the lifting surfaces. These can be a combination of Wing types constructed using the various methods available.

Info

Support for fuselages and control surfaces is in progress.

Horizontal tail

htail = WingSection(
    area       = 2.56 / 2,
    aspect     = 6.25,
    dihedral   = 0.0,
    sweep      = 15.0,
    taper      = 0.6,
    root_twist = 0.0,
    tip_twist  = 0.0,
    root_foil  = naca4(0,0,1,2),
    tip_foil   = naca4(0,0,0,9),
    position   = [5., 0., -0.1],
    angle      = 0.,
    axis       = [0., 1., 0.],
    symmetry   = true
);

Vertical tail

vtail = WingSection(
    area       = 0.512 / 2,
    aspect     = 1.25,
    dihedral   = 0.0,
    sweep      = 8.0,
    taper      = 0.6,
    root_twist = 0.0,
    tip_twist  = 0.0,
    root_foil  = naca4(0,0,0,9),
    tip_foil   = naca4(0,0,0,9),
    position   = [5., 0., 0.],
    angle      = 90.,
    axis       = [1., 0., 0.]
)
AbstractWing with 1 spanwise section(s).
Aspect Ratio: 1.2500000000000002
Span (m): 0.565685424949238
Projected Area (m²): 0.25599999999999995
Mean Aerodynamic Chord (m): 0.46197643037521097
Mean Aerodynamic Center (m): [5.151932479252752, 1.5875861030984043e-17, 0.2592724864350674]
Position (m): [5.0, 0.0, 0.0]

Meshing

The WingMesh type takes a Wing type with a vector of integers consisting of the spanwise panel distribution corresponding to the number of sections, and an integer for the chordwise distribution.

# Wing meshes
wing_mesh = WingMesh(wing, [20,8], 10)
WingMesh —
Symmetric AbstractWing with 2 spanwise section(s).
Aspect Ratio: 7.684336828804395
Span (m): 2.2
Projected Area (m²): 0.6298526610464908
Mean Aerodynamic Chord (m): 0.3037157904580297
Mean Aerodynamic Center (m): [0.10257127268682563, 0.0, 0.0]
Position (m): [0.0, 0.0, 0.0]
Spanwise panels: [4, 10, 10, 4]
Chordwise panels: 10
Spanwise spacing: AbstractSpacing[Cosine(), Sine(true), Sine(false), Cosine()]
Chordwise spacing: Cosine()
Tip

Optionally ,the type of spanwise spacing can be specified by the keyword span_spacing and providing types Sine(), Cosine(), Uniform() or a vector of the combination with length corresponding to the number of sections. For example, if you have two spanwise sections and the wing is symmetric, then the total number of spanwise sections is four. So the spanwise spacing argument would be, for example, spanwise_spacing = [Uniform(), Cosine(), Cosine(), Uniform()]

# Horizontal tail
htail_mesh = WingMesh(htail, [10], 5;
                      span_spacing = Cosine()) # Uniform(), Sine()

# Vertical tail
vtail_mesh = WingMesh(vtail, [10], 5); # Vertical tail

Inviscid Analysis

The inviscid 3D analysis uses a vortex lattice method. The WingMesh type allows you to generate horseshoe elements easily.

wing_horsies = make_horseshoes(wing_mesh)
10×28 Matrix{Horseshoe{Float64}}:
 Horseshoe{Float64}([0.147303, -1.1, 0.173152], [0.133232, -1.08536, 0.147788], [0.141581, -1.09268, 0.160365], [-5.14387e-6, 0.000129292, 7.75726e-5], 0.0)  …  Horseshoe{Float64}([0.133232, 1.08536, 0.147788], [0.147303, 1.1, 0.173152], [0.141581, 1.09268, 0.160365], [-5.14387e-6, -0.000129292, 7.75726e-5], 0.0)
 Horseshoe{Float64}([0.150914, -1.1, 0.172837], [0.137372, -1.08536, 0.147482], [0.147954, -1.09268, 0.159855], [7.70584e-7, 0.000386787, 0.000223864], 0.0)     Horseshoe{Float64}([0.137372, 1.08536, 0.147482], [0.150914, 1.1, 0.172837], [0.147954, 1.09268, 0.159855], [7.70584e-7, -0.000386787, 0.000223864], 0.0)
 Horseshoe{Float64}([0.159006, -1.1, 0.172132], [0.146649, -1.08536, 0.146798], [0.158763, -1.09268, 0.15899], [7.86514e-6, 0.000604229, 0.000348141], 0.0)      Horseshoe{Float64}([0.146649, 1.08536, 0.146798], [0.159006, 1.1, 0.172132], [0.158763, 1.09268, 0.15899], [7.86514e-6, -0.000604229, 0.000348141], 0.0)
 Horseshoe{Float64}([0.170787, -1.1, 0.171105], [0.160155, -1.08536, 0.145801], [0.17295, -1.09268, 0.157854], [2.21062e-5, 0.000763618, 0.000437704], 0.0)      Horseshoe{Float64}([0.160155, 1.08536, 0.145801], [0.170787, 1.1, 0.171105], [0.17295, 1.09268, 0.157854], [2.21062e-5, -0.000763618, 0.000437704], 0.0)
 Horseshoe{Float64}([0.185103, -1.1, 0.169857], [0.176568, -1.08536, 0.14459], [0.189127, -1.09268, 0.15656], [3.92157e-5, 0.000846632, 0.000484017], 0.0)       Horseshoe{Float64}([0.176568, 1.08536, 0.14459], [0.185103, 1.1, 0.169857], [0.189127, 1.09268, 0.15656], [3.92157e-5, -0.000846632, 0.000484017], 0.0)
 Horseshoe{Float64}([0.200554, -1.1, 0.168511], [0.194282, -1.08536, 0.143282], [0.205709, -1.09268, 0.155233], [4.81742e-5, 0.00084226, 0.000483298], 0.0)   …  Horseshoe{Float64}([0.194282, 1.08536, 0.143282], [0.200554, 1.1, 0.168511], [0.205709, 1.09268, 0.155233], [4.81742e-5, -0.00084226, 0.000483298], 0.0)
 Horseshoe{Float64}([0.215627, -1.1, 0.167197], [0.211562, -1.08536, 0.142007], [0.221074, -1.09268, 0.154003], [5.06927e-5, 0.000754789, 0.000435409], 0.0)     Horseshoe{Float64}([0.211562, 1.08536, 0.142007], [0.215627, 1.1, 0.167197], [0.221074, 1.09268, 0.154003], [5.06927e-5, -0.000754789, 0.000435409], 0.0)
 Horseshoe{Float64}([0.228846, -1.1, 0.166045], [0.226717, -1.08536, 0.140888], [0.233717, -1.09268, 0.152992], [4.52099e-5, 0.000594927, 0.000345142], 0.0)     Horseshoe{Float64}([0.226717, 1.08536, 0.140888], [0.228846, 1.1, 0.166045], [0.233717, 1.09268, 0.152992], [4.52099e-5, -0.000594927, 0.000345142], 0.0)
 Horseshoe{Float64}([0.238918, -1.1, 0.165167], [0.238264, -1.08536, 0.140036], [0.242402, -1.09268, 0.152297], [3.14006e-5, 0.0003798, 0.000221404], 0.0)       Horseshoe{Float64}([0.238264, 1.08536, 0.140036], [0.238918, 1.1, 0.165167], [0.242402, 1.09268, 0.152297], [3.14006e-5, -0.0003798, 0.000221404], 0.0)
 Horseshoe{Float64}([0.244856, -1.1, 0.164649], [0.245072, -1.08536, 0.139534], [0.246277, -1.09268, 0.151986], [1.12597e-5, 0.000130445, 7.62554e-5], 0.0)      Horseshoe{Float64}([0.245072, 1.08536, 0.139534], [0.244856, 1.1, 0.164649], [0.246277, 1.09268, 0.151986], [1.12597e-5, -0.000130445, 7.62554e-5], 0.0)
Note

You can also generate vortex ring elements by calling make_vortex rings, but the aerodynamic analysis is slightly slower.

For multiple lifting surfaces, it is most convenient to define a single vector consisting of all the components' horseshoes using ComponentArrays.jl.

aircraft = ComponentVector(
    wing = wing_horsies,
    htail = make_horseshoes(htail_mesh),
    vtail = make_horseshoes(vtail_mesh)
);

To define boundary conditions, use the following Freestream type, which takes named arguments for angles of attack and sideslip (in degrees), and a quasi-steady rotation vector.

# Define freestream conditions
fs = Freestream(
    alpha = 1.0,
    beta  = 0.0,
    omega = [0.,0.,0.]
)
Freestream: 
    alpha = 0.017453292519943295
    beta = 0.0
    omega = [0.0, 0.0, 0.0]

To define reference values, use the following References type.

# Define reference values
refs = References(
    speed     = 150.0,
    density   = 1.225,
    viscosity = 1.5e-5,
    area      = projected_area(wing),
    span      = span(wing),
    chord     = mean_aerodynamic_chord(wing),
    location  = mean_aerodynamic_center(wing)
)
References: 
    Mach = 0.45454545454545453
    Reynolds = 3.720518433110864e6
    speed = 150.0
    density = 1.225
    viscosity = 1.5e-5
    sound_speed = 330.0
    area = 0.6298526610464908
    span = 2.2
    chord = 0.3037157904580297
    location = [0.10257127268682563, 0.0, 0.0]

The vortex lattice analysis can be executed with the vortex elements, freestream conditions, and reference values defined. An optional named argument is provided for enabling compressibility corrections at $M > 0.3$ via the Prandtl-Glauert transformation.

# Solve system
system = solve_case(
    aircraft, fs, refs;
    compressible = true, # Compressibility corrections via Prandtl-Glauert transformation
    print = true, # Prints the results for only the aircraft
    print_components = true, # Prints the results for all components
)
VortexLatticeSystem -
380 Horseshoe{Float64} Elements
Freestream: 
    alpha = 0.017453292519943295
    beta = 0.0
    omega = [0.0, 0.0, 0.0]
References: 
    Mach = 0.45454545454545453
    Reynolds = 3.720518433110864e6
    speed = 150.0
    density = 1.225
    viscosity = 1.5e-5
    sound_speed = 330.0
    area = 0.6298526610464908
    span = 2.2
    chord = 0.3037157904580297
    location = [0.10257127268682563, 0.0, 0.0]

If needed, you can access the relevant component influence matrix values and boundary conditions with the following attributes, e.g.

system.influence_matrix[:wing] # Or :htail, :vtail
system.influence_matrix[:wing, :htail];

system.boundary_vector[:wing]
10×28 Matrix{Float64}:
  3.37518e-6   1.38378e-5    2.18839e-5   …   1.38378e-5    3.37518e-6
 -4.1663e-6   -1.45893e-6    1.07025e-5      -1.45893e-6   -4.1663e-6
 -1.24165e-5  -1.99414e-5   -5.74997e-6      -1.99414e-5   -1.24165e-5
 -2.64917e-5  -5.74585e-5   -4.82687e-5      -5.74585e-5   -2.64917e-5
 -4.24492e-5  -0.000102692  -0.00010299      -0.000102692  -4.24492e-5
 -5.04163e-5  -0.000126442  -0.000133125  …  -0.000126442  -5.04163e-5
 -5.19149e-5  -0.000133243  -0.000144428     -0.000133243  -5.19149e-5
 -4.56287e-5  -0.000118944  -0.00013137      -0.000118944  -4.56287e-5
 -3.14068e-5  -8.26604e-5   -9.23305e-5      -8.26604e-5   -3.14068e-5
 -1.12132e-5  -2.96488e-5   -3.32944e-5      -2.96488e-5   -1.12132e-5

Forces and Moments

The analysis computes the aerodynamic force coefficients $C_{\mathbf{F}}$ by surface pressure integration for nearfield forces. The aerodynamic moment coefficients $C_{\mathbf M}$ are calculated by computing the cross product of the moment arm of each control point from the reference point $\mathbf r_i - \mathbf r_{ref}$, i.e. $(C_{\mathbf{M}})_i = (\mathbf r_i - \mathbf r_{ref}) \times (C_{\mathbf{F}})_i$

Note

You can specify the axis system for the nearfield forces, with choices of geometry, body, wind, or stability axes. The wind axes are used by default.

# Compute dynamics
ax = Wind() # Body(), Stability(), Geometry()

# Compute non-dimensional coefficients
CFs, CMs = surface_coefficients(system; axes = ax)
((wing = StaticArraysCore.SVector{3, Float64}[[-1.1592060275150953e-6, 0.00013838990119929137, 8.136771701832143e-5] [-8.87423172689494e-7, 0.00040743538611431483, 0.0002380965040040056] … [-8.874231726894868e-7, -0.0004074353861143148, 0.00023809650400400556] [-1.1592060275150956e-6, -0.00013838990119929118, 8.13677170183213e-5]; [-1.1053927971495145e-6, 0.00020423816057873827, 0.0001197105394717196] [-1.152953183539807e-6, 0.0006241495203285039, 0.0003645952540669708] … [-1.1529531835398021e-6, -0.0006241495203285039, 0.0003645952540669708] [-1.1053927971495198e-6, -0.00020423816057873838, 0.00011971053947171966]; … ; [1.695193268127678e-6, 4.873646478769646e-5, 2.834398929143275e-5] [1.624139471604554e-6, 0.00013547800920556634, 7.89251592662278e-5] … [1.6241394716045505e-6, -0.00013547800920556604, 7.892515926622761e-5] [1.695193268127688e-6, -4.8736464787696593e-5, 2.8343989291432827e-5]; [3.299627834323949e-7, 1.0349195944237935e-5, 6.031562183009549e-6] [3.2551372798722074e-7, 2.828426237018824e-5, 1.6489314075997124e-5] … [3.255137279872168e-7, -2.8284262370188207e-5, 1.6489314075997104e-5] [3.2996278343239905e-7, -1.0349195944238008e-5, 6.0315621830095915e-6]], htail = StaticArraysCore.SVector{3, Float64}[[-7.444571948338659e-6, -1.1070992152192054e-5, 0.001975175895651409] [6.68360781424572e-6, -2.312909032497562e-5, 0.005400390813406512] … [6.683607814250848e-6, 2.3129090324973417e-5, 0.00540039081340633] [-7.444571948337348e-6, 1.107099215219155e-5, 0.001975175895651375]; [2.2758771166897212e-6, -6.462126637748469e-6, 0.0016562040870070737] [1.5497604853735847e-5, -1.833088014120721e-5, 0.005216093071606202] … [1.54976048537404e-5, 1.8330880141205383e-5, 0.005216093071606031] [2.2758771166905844e-6, 6.4621266377481125e-6, 0.001656204087007039]; … ; [3.202855879696291e-6, -2.509594969521728e-7, 0.00027733638646044014] [1.0702420574784663e-5, -1.635084555807346e-6, 0.0012245719609139396] … [1.0702420574785219e-5, 1.6350845558072273e-6, 0.001224571960913927] [3.2028558796963585e-6, 2.509594969521537e-7, 0.0002773363864604369]; [6.36514301065219e-7, -3.0529700831994e-8, 5.11157298694758e-5] [2.2915192676810637e-6, -2.040420391581354e-7, 0.0002291917799048206] … [2.291519267681131e-6, 2.0404203915812205e-7, 0.00022919177990481804] [6.365143010652399e-7, 3.052970083199214e-8, 5.111572986947611e-5]], vtail = StaticArraysCore.SVector{3, Float64}[[3.531765761310995e-29, -4.266766903499692e-16, -5.228616140714792e-30] [3.25509930076097e-28, -2.107538010800396e-15, -4.8302263455874217e-29] … [-2.809782151688737e-30, -2.5566400158096715e-16, 4.337520769372576e-31] [-5.439376071008263e-31, -5.3492112424312996e-17, 8.421426903980709e-32]; [-4.6321160894273365e-28, -1.8724307347132153e-15, 4.4086325049881833e-29] [-1.617974037588952e-27, -1.429753798708776e-14, 1.5446694671724906e-28] … [-2.9866368779511995e-30, -2.4679952454273537e-16, 2.9864868555288186e-31] [-3.510171491278477e-31, -5.651564888084964e-17, 3.678695624960685e-32]; … ; [-1.057369479708156e-25, 2.7622234835198218e-14, -1.3586716284456743e-26] [-1.8362086500119672e-27, 1.890197218045582e-14, -2.3706978807583864e-28] … [4.096615798879909e-32, 1.0182675941074257e-17, 4.641294980604185e-33] [-1.1266726991351777e-33, -4.0386551392468027e-17, 2.322287002255959e-33]; [-1.270033315287286e-24, -4.938741666901814e-14, -2.6987972788937287e-25] [5.675483316395856e-27, -1.9433353226887615e-14, 1.2072285131488992e-27] … [3.4081334252657876e-33, -5.120530835831756e-18, 1.0365616406819953e-33] [-1.2093024186508735e-31, -4.9317597445476505e-17, -2.2689566846762384e-32]]), (wing = StaticArraysCore.SizedVector{3, Float64, Vector{Float64}}[[5.046439601658936e-5, -1.1457792138756693e-5, -1.971341031459235e-6] [0.0001372025098669018, -1.2664801161380094e-5, -2.4805364809174786e-6] … [-0.00013720250986690177, -1.2664801161380089e-5, 2.480536480917482e-6] [-5.046439601658928e-5, -1.1457792138756676e-5, 1.9713410314592316e-6]; [7.425573075422097e-5, -1.8065047228930806e-5, -3.5692210555022104e-6] [0.00021000304032203608, -2.5041487285910542e-5, -5.254012986445412e-6] … [-0.00021000304032203608, -2.5041487285910535e-5, 5.254012986445414e-6] [-7.425573075422102e-5, -1.8065047228930813e-5, 3.5692210555022104e-6]; … ; [1.7405127102710345e-5, -1.2102121601405733e-5, -3.913722740533952e-6] [4.491245703159879e-5, -3.497485676439708e-5, -9.212299297039867e-6] … [-4.491245703159869e-5, -3.4974856764397005e-5, 9.212299297039847e-6] [-1.740512710271039e-5, -1.2102121601405766e-5, 3.9137227405339656e-6]; [3.6993729382272335e-6, -2.717577008462002e-6, -8.461073928069136e-7] [9.374287576908095e-6, -7.738695682856828e-6, -2.0176021082505015e-6] … [-9.374287576908086e-6, -7.73869568285682e-6, 2.0176021082504977e-6] [-3.6993729382272598e-6, -2.71757700846202e-6, 8.461073928069204e-7]], htail = StaticArraysCore.SizedVector{3, Float64, Vector{Float64}}[[0.00121003448968918, -0.0342295952858522, 3.104736357990523e-5] [0.002708071499507066, -0.09245809827150885, 5.131516760432633e-5] … [-0.0027080714995069747, -0.09245809827150571, -5.1315167604318566e-5] [-0.001210034489689159, -0.034229595285851615, -3.104736357990323e-5]; [0.001014381958571387, -0.02896333165891721, 1.4207170062013081e-5] [0.0026153138030529234, -0.09020510948784716, 3.5993283778778626e-5] … [-0.0026153138030528375, -0.0902051094878442, -3.599328377877199e-5] [-0.001014381958571366, -0.0289633316589166, -1.4207170062011688e-5]; … ; [0.00016978891037332732, -0.0050305686899799786, -1.3323972666991528e-6] [0.0006137624267870493, -0.022059199883957024, -1.2978988876621612e-6] … [-0.0006137624267870431, -0.022059199883956795, 1.2978988876627353e-6] [-0.0001697889103733253, -0.005030568689979919, 1.3323972666992422e-6]; [3.129232469507141e-5, -0.0009396708804145746, -3.1218523423880137e-7] [0.00011486329423311801, -0.00419089963405393, -6.333557235899161e-7] … [-0.00011486329423311669, -0.004190899634053884, 6.333557235899826e-7] [-3.12923246950716e-5, -0.0009396708804145802, 3.121852342388189e-7]], vtail = StaticArraysCore.SizedVector{3, Float64, Vector{Float64}}[[1.5283331633297508e-17, 7.538532785212894e-29, 9.524999807814655e-16] [4.9680046429039635e-17, 7.261245864490791e-28, 4.708643649120535e-15] … [-5.168485422110365e-17, -1.1239747939219533e-29, 5.797890749575407e-16] [-1.1469041144467642e-17, -2.22926775495351e-30, 1.214056822746722e-16]; [6.81908996623334e-17, -6.016562509140591e-28, 4.24420365666465e-15] [3.4542811860601737e-16, -2.2543349460843615e-27, 3.242459888627923e-14] … [-4.980013043461017e-17, -9.317816920529054e-30, 5.649964685883181e-16] [-1.2096730897095959e-17, -1.1545698037562833e-30, 1.2944665458390947e-16]; … ; [-1.074911315105466e-15, 2.6696033433105676e-25, -6.656118479387956e-14] [-5.029510664250239e-16, 4.489228499057619e-27, -4.5518117924389565e-14] … [2.0387589316892223e-18, -2.0567036640971524e-32, -2.4224302550445862e-17] [-8.58313996273349e-18, -4.172633312679413e-32, 9.601462571785542e-17]; [1.968513891326863e-15, 5.1831235263290956e-24, 1.2167932748271943e-13] [5.350813610169349e-16, -2.265399564694757e-26, 4.782845388292172e-14] … [-1.022193655616286e-18, -1.318878361880758e-32, 1.2355232331168076e-17] [-1.0452911034921555e-17, 2.104742713491705e-31, 1.188683409516016e-16]]))

You can access the corresponding values of the components' by the name provided in the ComponentVector.

CFs.wing
10×28 reshape(view(::Vector{StaticArraysCore.SVector{3, Float64}}, 1:280), 10, 28) with eltype StaticArraysCore.SVector{3, Float64}:
 [-1.15921e-6, 0.00013839, 8.13677e-5]    …  [-1.15921e-6, -0.00013839, 8.13677e-5]
 [-1.10539e-6, 0.000204238, 0.000119711]     [-1.10539e-6, -0.000204238, 0.000119711]
 [1.2262e-6, 0.000203867, 0.000118253]       [1.2262e-6, -0.000203867, 0.000118253]
 [4.63886e-6, 0.000209934, 0.00012037]       [4.63886e-6, -0.000209934, 0.00012037]
 [6.94506e-6, 0.000209329, 0.00011958]       [6.94506e-6, -0.000209329, 0.00011958]
 [7.1121e-6, 0.000184144, 0.000105474]    …  [7.1121e-6, -0.000184144, 0.000105474]
 [5.93898e-6, 0.000149396, 8.60475e-5]       [5.93898e-6, -0.000149396, 8.60475e-5]
 [3.84756e-6, 0.000101061, 5.85346e-5]       [3.84756e-6, -0.000101061, 5.85346e-5]
 [1.69519e-6, 4.87365e-5, 2.8344e-5]         [1.69519e-6, -4.87365e-5, 2.8344e-5]
 [3.29963e-7, 1.03492e-5, 6.03156e-6]        [3.29963e-7, -1.03492e-5, 6.03156e-6]

Functions are available for directly retrieving the dimensionalized forces and moments.

Fs, Ms = surface_dynamics(system; axes = ax)
((wing = StaticArraysCore.SVector{3, Float64}[[-0.010062090296843686, 1.2012460675549868, 0.7062845572339066] [-0.007702972451113239, 3.5366031127354725, 2.066715032345629] … [-0.007702972451113177, -3.5366031127354716, 2.0667150323456287] [-0.010062090296843688, -1.201246067554985, 0.7062845572339055]; [-0.009594983009397949, 1.7728192961606764, 1.039106275379186] [-0.010007814629535178, 5.417715818593168, 3.164744041304385] … [-0.010007814629535136, -5.417715818593168, 3.164744041304385] [-0.009594983009397994, -1.7728192961606772, 1.0391062753791864]; … ; [0.014714543687343037, 0.42304016525341825, 0.246030276632126] [0.014097785579138739, 1.1759703879259522, 0.6850827725006905] … [0.014097785579138706, -1.1759703879259495, 0.6850827725006889] [0.014714543687343123, -0.4230401652534195, 0.2460302766321267]; [0.002864128759416238, 0.08983264546499652, 0.05235490661359539] [0.002825510259716948, 0.24551183757949013, 0.14312983475539948] … [0.002825510259716914, -0.24551183757948986, 0.1431298347553993] [0.002864128759416274, -0.08983264546499714, 0.05235490661359576]], htail = StaticArraysCore.SVector{3, Float64}[[-0.06462005319805642, -0.09609794986124985, 17.14483684733495] [0.05801476505413938, -0.20076413494211903, 46.87624004097393] … [0.0580147650541839, 0.2007641349420999, 46.87624004097235] [-0.06462005319804504, 0.0960979498612455, 17.144836847334656]; [0.019754970651542803, -0.05609227367291036, 14.376111474497764] [0.1345216430226528, -0.15911491729975824, 45.27650671015759] … [0.13452164302269232, 0.1591149172997424, 45.276506710156106] [0.019754970651550297, 0.05609227367290726, 14.376111474497462]; … ; [0.027801291836244465, -0.0021783678304332726, 2.4073233721422707] [0.09289869070912732, -0.01419279062823353, 10.629476860219777] … [0.09289869070913215, 0.0141927906282325, 10.629476860219668] [0.02780129183624505, 0.0021783678304331065, 2.4073233721422422]; [0.005525044056473539, -0.00026500259592822655, 0.44369255967230353] [0.019890746977731466, -0.001771116931441713, 1.9894206292561725] … [0.01989074697773205, 0.0017711169314415971, 1.9894206292561503] [0.005525044056473721, 0.0002650025959282104, 0.44369255967230625]], vtail = StaticArraysCore.SVector{3, Float64}[[3.0656281242593306e-25, -3.703620654098001e-12, -4.538520891595474e-26] [2.8254772932521786e-24, -1.829376078570102e-11, -4.192712295300848e-25] … [-2.438935017044125e-26, -2.2192036691480875e-12, 3.76503612040556e-27] [-4.721463819705962e-27, -4.643199333048174e-13, 7.30993075646505e-28]; [-4.020749482919107e-24, -1.6252992720937516e-11, 3.8267622232678228e-25] [-1.4044268644002363e-23, -1.2410487422759396e-10, 1.3407973465066049e-24] … [-2.592447695794693e-26, -2.1422586168658914e-12, 2.5923174739769386e-27] [-3.046883958873355e-27, -4.90564704397568e-13, 3.1931655524864157e-28]; … ; [-9.178133075064212e-22, 2.397653346473531e-10, -1.1793483026037838e-22] [-1.5938579339405014e-23, 1.6407208585334678e-10, -2.0578029769100916e-24] … [3.5559268241701146e-28, 8.838722569638527e-14, 4.0287169045554783e-29] [-9.779695899259917e-30, -3.505616061712891e-13, 2.0157815743915813e-29]; [-1.1024088552933275e-20, -4.2869052977300213e-10, -2.342598405161511e-21] [4.926408615233093e-23, -1.6868455675525328e-10, 1.0478933010957265e-23] … [2.9583133157292823e-29, -4.444701150179333e-14, 8.997517765797551e-30] [-1.0496934836286554e-27, -4.2808448795208666e-13, -1.9694900215261427e-28]]), (wing = StaticArraysCore.SVector{3, Float64}[[0.9636855350750623, -0.03020618589432311, -0.0376454091730872] [2.620066513651755, -0.03338822467386794, -0.047369181335299215] … [-2.6200665136517545, -0.03338822467386793, 0.04736918133529927] [-0.9636855350750609, -0.030206185894323064, 0.03764540917308713]; [1.4180130799692376, -0.04762489737800592, -0.06815907796740718] [4.010290586131321, -0.06601689146288019, -0.1003324465524074] … [-4.010290586131321, -0.06601689146288017, 0.10033244655240744] [-1.4180130799692385, -0.047624897378005945, 0.06815907796740718]; … ; [0.33237431831168673, -0.031904832133517164, -0.07473780112432236] [0.8576637907605976, -0.09220424076984002, -0.17592124900139625] … [-0.8576637907605956, -0.09220424076983981, 0.17592124900139586] [-0.3323743183116876, -0.03190483213351725, 0.07473780112432263]; [0.07064450327010295, -0.00716435027845157, -0.016157558990700156] [0.17901463313918292, -0.020401529155453364, -0.03852882667254731] … [-0.17901463313918273, -0.020401529155453346, 0.03852882667254724] [-0.07064450327010345, -0.0071643502784516176, 0.016157558990700285]], htail = StaticArraysCore.SVector{3, Float64}[[23.107236521211178, -90.23950738244872, 0.5928911776598951] [51.714268633396244, -243.74735289340333, 0.9799321631430065] … [-51.7142686333945, -243.7473528933951, -0.979932163142858] [-23.107236521210776, -90.23950738244717, -0.5928911776598567]; [19.370988215037876, -76.35605268565476, 0.27130502619337216] [49.942935626451, -237.80779689579433, 0.6873401779357136] … [-49.94293562644936, -237.80779689578654, -0.6873401779355868] [-19.37098821503747, -76.35605268565315, -0.27130502619334557]; … ; [3.2423476719932065, -13.262091960082072, -0.025443918371071465] [11.720619275276945, -58.15468497817771, -0.024785125410375836] … [-11.720619275276826, -58.154684978177116, 0.0247851254103868] [-3.242347671993168, -13.262091960081916, 0.025443918371073175]; [0.597569039716619, -2.477255037405914, -0.005961597051534192] [2.1934691366788264, -11.04847180657794, -0.012094779637904653] … [-2.1934691366788006, -11.048471806577819, 0.012094779637905922] [-0.5975690397166227, -2.477255037405929, 0.005961597051534527]], vtail = StaticArraysCore.SVector{3, Float64}[[2.9185577922942417e-13, 1.9873839560271451e-25, 1.818926859516212e-11] [9.487073244626678e-13, 1.914282784597576e-24, 8.991788533422331e-11] … [-9.86991826454582e-13, -2.963135580972768e-26, 1.107185241546254e-11] [-2.1901677072427177e-13, -5.8770202320717794e-27, 2.3184048382459516e-12]; [1.3021969708454063e-12, -1.586146819518376e-24, 8.104877883599087e-11] [6.596414652409228e-12, -5.943104886600098e-24, 6.191913387221813e-10] … [-9.510004901061078e-13, -2.4564567642899413e-26, 1.0789367695356252e-11] [-2.3100335102385375e-13, -3.0437932280397274e-27, 2.471958022550461e-12]; … ; [-2.0526877712261182e-11, 7.037877269511691e-22, -1.2710753728675086e-9] [-9.604527267205168e-12, 1.1834956414154501e-23, -8.692297003446783e-10] … [3.893284467980395e-14, -5.422089391647437e-29, -4.625956477805155e-13] [-1.6390660506254486e-13, -1.1000316290020981e-28, 1.8335284530013256e-12]; [3.759142112871823e-11, 1.366427238804186e-20, 2.3236304616468746e-9] [1.0218098469479072e-11, -5.972274548830649e-23, 9.13348673722846e-10] … [-1.9520163079710448e-14, -3.476960001413802e-29, 2.3593978368678637e-13] [-1.996123992144672e-13, 5.548734773115831e-28, 2.2699508920263724e-12]]))

A Trefftz plane integration is performed to compute farfield forces.

Note

The farfield forces are usually more accurate compared to nearfield forces.

To obtain the nearfield coefficients of the components (in wind axes by definition):

nfs = nearfield_coefficients(system)
(wing = 6-element LabelledArrays.SLArray{Tuple{6}, Float64, 1, 6, (:CX, :CY, :CZ, :Cl, :Cm, :Cn)} with indices SOneTo(6):
 :CX => 0.003957539725248549
 :CY => -1.0477797557535695e-17
 :CZ => 0.3450556800095906
 :Cl => 2.540040050578833e-17
 :Cm => -0.05858161327748336
 :Cn => 3.0535537748517527e-19, htail = 6-element LabelledArrays.SLArray{Tuple{6}, Float64, 1, 6, (:CX, :CY, :CZ, :Cl, :Cm, :Cn)} with indices SOneTo(6):
 :CX => 0.0003325950091132996
 :CY => -1.2038727078092661e-17
 :CZ => 0.053102126952264656
 :Cl => 2.419600435808744e-16
 :Cm => -0.9209847143167378
 :Cn => 4.1892873142658407e-17, vtail = 6-element LabelledArrays.SLArray{Tuple{6}, Float64, 1, 6, (:CX, :CY, :CZ, :Cl, :Cm, :Cn)} with indices SOneTo(6):
 :CX => -1.3837498955910852e-24
 :CY => -1.0340285225910211e-14
 :CZ => -2.948680861186249e-25
 :Cl => -1.6833873814293458e-16
 :Cm => 5.636762354961024e-24
 :Cn => 1.8795165493456688e-14)

Similarly for the farfield coefficients of the components.

ffs = farfield_coefficients(system)
(wing = 3-element LabelledArrays.SLArray{Tuple{3}, Float64, 1, 3, (:CDi, :CY, :CL)} with indices SOneTo(3):
 :CDi => 0.004528682767972114
  :CY => -4.112899017506603e-18
  :CL => 0.34503920152404866, htail = 3-element LabelledArrays.SLArray{Tuple{3}, Float64, 1, 3, (:CDi, :CY, :CL)} with indices SOneTo(3):
 :CDi => 0.00035475639005410026
  :CY => 6.6662473193905846e-18
  :CL => 0.053335728251027344, vtail = 3-element LabelledArrays.SLArray{Tuple{3}, Float64, 1, 3, (:CDi, :CY, :CL)} with indices SOneTo(3):
 :CDi => 1.951432706334011e-25
  :CY => -1.0338756055059793e-14
  :CL => 6.330662254997145e-31)

You can access the values corresponding to the components by the name used in the ComponentArray construction.

@show (nfs.wing, ffs.wing)
(6-element LabelledArrays.SLArray{Tuple{6}, Float64, 1, 6, (:CX, :CY, :CZ, :Cl, :Cm, :Cn)} with indices SOneTo(6):
 :CX => 0.003957539725248549
 :CY => -1.0477797557535695e-17
 :CZ => 0.3450556800095906
 :Cl => 2.540040050578833e-17
 :Cm => -0.05858161327748336
 :Cn => 3.0535537748517527e-19, 3-element LabelledArrays.SLArray{Tuple{3}, Float64, 1, 3, (:CDi, :CY, :CL)} with indices SOneTo(3):
 :CDi => 0.004528682767972114
  :CY => -4.112899017506603e-18
  :CL => 0.34503920152404866)

To obtain the total nearfield and farfield force coefficients:

nf, ff = nearfield(system), farfield(system)
(6-element LabelledArrays.SLArray{Tuple{6}, Float64, 1, 6, (:CX, :CY, :CZ, :Cl, :Cm, :Cn)} with indices SOneTo(6):
 :CX => 0.0042901347343618475
 :CY => -1.0362800056479945e-14
 :CZ => 0.3981578069618553
 :Cl => 9.79074341016154e-17
 :Cm => -0.9795663275942211
 :Cn => 1.883735652219678e-14, 3-element LabelledArrays.SLArray{Tuple{3}, Float64, 1, 3, (:CDi, :CY, :CL)} with indices SOneTo(3):
 :CDi => 0.004883439158026214
  :CY => -1.0336202706757908e-14
  :CL => 0.398374929775076)

You can also print the relevant information as a pretty table, if necessary.

print_coefficients(nfs.wing, ffs.wing, :wing)
print_coefficients(nf, ff, :aircraft)
────────────────────────────────────
 wing  Nearfield          Farfield
────────────────────────────────────
  CX   0.00395754  CDi   0.00452868
  CY      -0.0     CYff     -0.0
  CZ    0.345056    CL    0.345039
  Cl      0.0                —
  Cm   -0.0585816            —
  Cn      0.0                —
────────────────────────────────────
────────────────────────────────────────
 aircraft  Nearfield          Farfield
────────────────────────────────────────
    CX     0.00429013  CDi   0.00488344
    CY        -0.0     CYff     -0.0
    CZ      0.398158    CL    0.398375
    Cl        0.0                —
    Cm     -0.979566             —
    Cn        0.0                —
────────────────────────────────────────

Aerodynamic Stability Analyses

The derivatives of the aerodynamic coefficients with respect to the freestream values are obtained by automatic differentiation enabled by ForwardDiff. The following function evaluates the derivatives of the aerodynamic coefficients with respect to the freestream values. You can also optionally provide the axes for the reference frame of the coefficients.

dv_data = freestream_derivatives(
    system,
    print = true,  # Prints the results for only the aircraft
    print_components = true,  # Prints the results for all components
);
──────────────────────────────────────────────────────────────────────────────────────────────
 wing    Values                             Freestream   Derivatives
                      ∂/∂M     ∂/∂α, 1/rad  ∂/∂β, 1/rad     ∂/∂p̄         ∂/∂q̄        ∂/∂r̄
──────────────────────────────────────────────────────────────────────────────────────────────
  CX   0.00395754  0.00347335    0.11771       -0.0          0.0      -0.0492756      0.0
  CY      -0.0        0.0         -0.0       -0.142373    -0.264392      0.0       0.0918817
  CZ    0.345056    0.151508      5.472        -0.0          0.0       10.9144        0.0
  Cℓ      0.0         0.0         -0.0      -0.0418825    -0.555361      -0.0      0.109052
  Cm   -0.0585816  -0.0319284   0.0299961       0.0         -0.0      -0.906578      -0.0
  Cn      -0.0        -0.0         0.0      -0.0080398   -0.0253157      0.0      -0.00233701
──────────────────────────────────────────────────────────────────────────────────────────────
─────────────────────────────────────────────────────────────────────────────────────────────
 htail   Values                              Freestream   Derivatives
                      ∂/∂M      ∂/∂α, 1/rad  ∂/∂β, 1/rad     ∂/∂p̄        ∂/∂q̄       ∂/∂r̄
─────────────────────────────────────────────────────────────────────────────────────────────
  CX    0.0003326  -0.00052487    0.12502       -0.0          0.0      3.72439       0.0
  CY      -0.0         0.0         -0.0      0.00230357    0.014003      -0.0    0.00244072
  CZ    0.0531021  -0.0343939     8.53404       -0.0         -0.0      403.809       0.0
  Cℓ       0.0        -0.0          0.0       -0.203665    -1.42212      0.0      -0.52356
  Cm    -0.920985   0.577482     -145.723        0.0          0.0      -6896.25     -0.0
  Cn       0.0        -0.0          0.0      -0.0069241    -0.038453     0.0     -0.00135864
─────────────────────────────────────────────────────────────────────────────────────────────
────────────────────────────────────────────────────────────────────────────
 vtail  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     -2.0e-8     -1.02913     0.332084    0.0   4.76381
  CZ     -0.0   -0.0     -0.0         -0.0          0.0      0.0     0.0
  Cℓ     -0.0   0.0      -0.0      -0.0875382   0.00660004   0.0   0.441621
  Cm     0.0    0.0       0.0          0.0         -0.0      -0.0    -0.0
  Cn     0.0    -0.0    3.0e-8       2.30749     -0.611054   -0.0  -10.9476
────────────────────────────────────────────────────────────────────────────
──────────────────────────────────────────────────────────────────────────────────────────────
 aircraft    Values                             Freestream   Derivatives
                          ∂/∂M     ∂/∂α, 1/rad  ∂/∂β, 1/rad     ∂/∂p̄        ∂/∂q̄      ∂/∂r̄
──────────────────────────────────────────────────────────────────────────────────────────────
    CX     0.00429013  0.00294848    0.24273       -0.0          0.0      3.67512      0.0
    CY        -0.0        0.0        -2.0e-8      -1.1692     0.0816947     -0.0     4.85813
    CZ      0.398158    0.117114     14.006        -0.0          0.0      414.723      0.0
    Cℓ        0.0         -0.0         0.0       -0.333086    -1.97088      0.0     0.0271125
    Cm     -0.979566    0.545553    -145.693        0.0         -0.0      -6897.15    -0.0
    Cn        0.0         -0.0       3.0e-8       2.29253     -0.674823     0.0     -10.9513
──────────────────────────────────────────────────────────────────────────────────────────────
Note

For efficiency, instead of calling solve_case to compute the forces and then computing the derivatives, you can directly call:

dvs = freestream_derivatives(
          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
      )

This page was generated using Literate.jl.