Aerodynamics API

Doublet-Source Panel API

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

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

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

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

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

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

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

source
AeroFuse.DoubletSource.source_strengthsMethod
source_strengths(panels, freestream)

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

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

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

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

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

source

Vortex Lattice API

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

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

The finite core setup is not implemented for now.

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

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

Arguments

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

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

source
AeroFuse.VortexLattice.VortexLatticeSystemType
VortexLatticeSystem

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

Arguments

The accessible fields are:

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

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

    r1 —front leg→ r4
    |               |
left leg       right leg
    ↓               ↓
    r2 —back leg-→ r3
source
AeroFuse.VortexLattice.boundary_conditionMethod
boundary_condition(vortices, U, Ω)

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

Compute the velocity of Freestream in the geometry axis system.

source