Library
Public API
Generating Lifting Surfaces
VortexLattice.AbstractSpacing
— TypeAbstractSpacing
Spacing discretization scheme supertype
VortexLattice.Uniform
— TypeUniform()
Uniform discretization scheme.
VortexLattice.Sine
— TypeSine()
Sine-spaced discretization scheme. Using sine-spacing on the right half of a wing effectively results in cosine spacing once symmetry is applied.
VortexLattice.Cosine
— TypeCosine()
Cosine-spaced discretization scheme. This is typically one of the most accurate spacing schemes for spanwise spacing.
VortexLattice.SurfacePanel
— TypeSurfacePanel{TF}
Lifting surface panel with attached vortex ring
Fields
rtl
: position of the left side of the top bound vortexrtc
: position of the center of the top bound vortexrtr
: position of the right side of the top bound vortexrbl
: position of the left side of the bottom bound vortexrbc
: position of the center of the bottom bound vortexrbr
: position of the right side of the bottom bound vortexrcp
: position of the panel control pointncp
: normal vector at the panel control pointcore_size
: finite core size (for use when the finite core smoothing model is enabled)chord
: panel chord length (for determining unsteady forces)
VortexLattice.SurfacePanel
— MethodSurfacePanel(rtl, rtr, rbl, rbr, rcp, ncp, core_size, chord; kwargs...)
Construct and return a vortex ring panel.
Arguments
rtl
: position of the left side of the top bound vortexrtr
: position of the right side of the top bound vortexrbl
: position of the left side of the bottom bound vortexrbr
: position of the right side of the bottom bound vortexrcp
: position of the panel control pointncp
: normal vector at the panel control pointcore_size
: finite core size (for use when the finite core smoothing model is enabled)chord
: panel chord length (for determining unsteady forces)
Keyword Arguments
rtc
: position of the center of the top bound vortex, defaults to(rtl+rtr)/2
rbc
: position of the center of the bottom bound vortex, defaults to(rbl+rbr)/2
VortexLattice.WakePanel
— TypeWakePanel{TF}
SurfacePanel used for modeling wakes.
Fields
rtl
: position of the left side of the top bound vortexrtr
: position of the right side of the top bound vortexrbl
: position of the left side of the bottom bound vortexrbr
: position of the right side of the bottom bound vortexcore_size
: finite core size (for use when the finite core smoothing model is enabled)gamma
: circulation strength of the panel
VortexLattice.WakePanel
— MethodWakePanel(rtl, rtr, rbl, rbr, core_size, gamma)
Construct and return a wake panel.
Arguments
rtl
: position of the left side of the top bound vortexrtr
: position of the right side of the top bound vortexrbl
: position of the left side of the bottom bound vortexrbr
: position of the right side of the bottom bound vortexcore_size
: finite core sizegamma
: circulation strength of the panel
VortexLattice.grid_to_surface_panels
— Functiongrid_to_surface_panels(xyz; mirror = false, fcore = (c, Δs) -> 1e-3)
Construct a set of panels with associated vortex rings given a potentially curved lifting surface defined by a grid with dimensions (3, i, j) where i
corresponds to the chordwise direction (ordered from leading edge to trailing edge) and j
corresponds to the spanwise direction (ordered from left to right). The leading edge of each ring vortex will be placed at the 1/4 chord and the control point will be placed at the 3/4 chord of each panel.
Return a grid with dimensions (3, i, j) containing the panel corners and a matrix with dimensions (i, j) containing the generated panels.
Keyword Arguments
mirror
: mirror the geometry across the X-Z plane? defaults tofalse
.fcore
: function for setting the finite core size based on the chord length (in the x-direction) and/or the panel width (in the y/z directions). Defaults to(c, Δs) -> 1e-3
grid_to_surface_panels(xyz, ns, nc;
mirror = false,
fcore = (c, Δs) -> 1e-3,
spacing_s = Cosine(),
spacing_c = Uniform(),
interp_s = (x, y, xpt) -> FLOWMath.linear(x, y, xpt),
interp_c = (x, y, xpt) -> FLOWMath.linear(x, y, xpt))
Discretize a potentially curved lifting surface defined by a grid with dimensions (3, i, j) where i
corresponds to the chordwise direction (ordered from leading edge to trailing edge) and j
corresponds to the spanwise direction (ordered from left to right) into ns
spanwise and nc
chordwise panels with associated vortex rings according to the spanwise discretization scheme spacing_s
and chordwise discretization scheme spacing_c
. The bound vortex will be placed at the 1/4 chord and the control point will be placed at the 3/4 chord of each panel.
Return a grid with dimensions (3, i, j) containing the interpolated panel corners and a matrix with dimensions (i, j) containing the generated panels.
Arguments
xyz
: grid of dimensions (3, i, j) where wherei
corresponds to the chordwise direction andj
corresponds to the spanwise direction.ns
: number of spanwise panelsnc
: number of chordwise panels
Keyword Arguments
mirror
: mirror the geometry across the X-Z plane? defaults tofalse
.fcore
: function for setting the finite core size based on the chord length (in the x-direction) and/or the panel width (in the y/z directions). Defaults to(c, Δs) -> 1e-3
spacing_s
: spanwise discretization scheme, defaults toCosine()
spacing_c
: chordwise discretization scheme, defaults toUniform()
interp_s
: spanwise interpolation function, defaults to linear interpolationinterp_c
: chordwise interpolation function, defaults to linear interpolation
VortexLattice.wing_to_surface_panels
— Functionwing_to_surface_panels(xle, yle, zle, chord, theta, phi, ns, nc;
fc = fill(x -> 0, length(xle)),
mirror = false,
fcore = (c, Δs) -> 1e-3,
spacing_s = Cosine(),
spacing_c = Uniform(),
interp_s = (x, y, xpt) -> FLOWMath.linear(x, y, xpt))
Discretize a wing into ns
spanwise and nc
chordwise panels with associated vortex rings according to the spanwise discretization scheme spacing_s
and chordwise discretization scheme spacing_c
.
Return a grid with dimensions (3, i, j) containing the panel corners and a matrix with dimensions (i, j) containing the generated panels.
Arguments
xle
: leading edge x-coordinate of each airfoil sectionyle
: leading edge y-coordinate of each airfoil sectionzle
: leading edge z-coordinate of each airfoil sectionchord
: chord length of each airfoil sectiontheta
: twist of each airfoil sectionphi
: dihedral angle of each airfoil section, defined by a right hand rotation about the x-axisns
: number of spanwise panelsnc
: number of chordwise panelsfc
: (optional) camber line function y=f(x) of each airfoil section- 'reference_line': 2D array, each row is the x, y coordinate of the reference point of the airfoil. This allows xle, yle, and zle to be defined about points that are not the leading edge
mirror
: mirror the geometry across the X-Z plane?, defaults tofalse
fcore
: function for setting the finite core size based on the chord length (in the x-direction) and/or the panel width (in the y/z directions). Defaults to(c, Δs) -> 1e-3
spacing_s
: spanwise discretization scheme, defaults toCosine()
spacing_c
: chordwise discretization scheme, defaults toUniform()
interp_s
: interpolation function between spanwise stations, defaults to linear interpolation
VortexLattice.lifting_line_geometry
— Functionlifting_line_geometry(grids, xc=0.25)
Construct a lifting line representation of the surfaces in grids
at the normalized chordwise location xc
. Return the lifting line coordinates and chord lengths.
Arguments
grids
: Vector with length equal to the number of surfaces. Each element of the vector is a grid with shape (3, nc, ns) which defines the discretization of a surface into panels.nc
is the number of chordwise panels andns
is the number of spanwise panels.xc
: Normalized chordwise location of the lifting line from the leading edge. Defaults to the quarter chord
Return Arguments:
r
: Vector with length equal to the number of surfaces, with each element being a matrix with size (3, ns+1) which contains the x, y, and z coordinates of the resulting lifting line coordinatesc
: Vector with length equal to the number of surfaces, with each element being a vector of lengthns+1
which contains the chord lengths at each lifting line coordinate.
VortexLattice.lifting_line_geometry!
— Functionlifting_line_geometry!(r, c, grids, xc=0.25)
In-place version of lifting_line_geometry
VortexLattice.read_degengeom
— Function`read_degengeom(filename::String)`
Read all geometry components from a DegenGeom file written out by OpenVSP
Arguments
filename::String
: DegenGeom filename
Returns
comp
: Vector ofvsp.VSPComponent
objects
VortexLattice.import_vsp
— Function`import_vsp(comp::vsp.VSPComponent; geomType::String="", optargs...)
Imports properties from OpenVSP component to VortexLattice objects. Importing prop and duct geometries are under development.
Arguments
comp::VSPGeom.VSPComponent
: SingleVSPGeom.VSPComponent
objectgeomType::String
: Geometry type may be one of -wing
,fuselage
,prop
,duct
optargs
: Optional arguments that are passed into gridtosurface_panels() called inside
Returns
grid
: Array with dimensions (3, i, j) containing the panel cornerssurface
: Array with dimensions (i, j) containing generated panels
VortexLattice.set_normal
— Functionset_normal(panel::SurfacePanel, ncp)
Return a copy of panel
with the new normal vector ncp
VortexLattice.translate
— Functiontranslate(panel::SurfacePanel, r)
Return a copy of panel
translated the distance specified by vector r
translate(surface, r)
Return a copy of the panels in surface
translated the distance specified by vector r
translate(surfaces, r)
Return a copy of the surfaces in surfaces
translated the distance specified by vector r
translate(grid, r)
Return a copy of the grid points in grid
translated the distance specified by vector r
VortexLattice.translate!
— Functiontranslate!(surface, r)
Translate the panels in surface
the distance specified by vector r
translate!(surfaces, r)
Translate the surfaces in surfaces
the distance specified by vector r
translate!(grid, r)
Translate the grid points in grid
the distance specified by vector r
VortexLattice.rotate
— Functionrotate(panel::SurfacePanel, R, r = [0,0,0])
Return a copy of panel
rotated about point r
using the rotation matrix R
rotate(surface, R, r = [0,0,0])
Return a copy of the panels in surface
rotated about point r
using the rotation matrix R
rotate(surfaces, R, r = [0,0,0])
Return a copy of the surfaces in surfaces
rotated about point r
using the rotation matrix R
rotate(grid, R, r = [0,0,0])
Return a copy of the grid points in grid
rotated about point r
using the rotation matrix R
VortexLattice.rotate!
— Functionrotate!(surface, R, r = [0,0,0])
Rotate the panels in surface
about point r
using the rotation matrix R
rotate!(surfaces, R, r = [0,0,0])
Rotate the surfaces in surfaces
about point r
using the rotation matrix R
rotate!(grid, R, r = [0,0,0])
Rotate the grid points in grid
about point r
using the rotation matrix R
VortexLattice.reflect
— Methodreflect(surface)
Reflects the panels in surface
across the X-Z plane
Reference Parameters and Frames
VortexLattice.Reference
— TypeReference(S, c, b, r)
Reference quantities.
Arguments
S
: reference areac
: reference chordb
: reference spanr
: reference location for all rotations/momentsV
: reference velocity (magnitude)
VortexLattice.AbstractFrame
— TypeAbstractFrame
Supertype for the different possible reference frames used by this package.
VortexLattice.Body
— TypeBody <: AbstractFrame
Reference frame aligned with the global X-Y-Z axes
VortexLattice.Stability
— TypeStability <: AbstractFrame
Reference frame rotated from the body frame about the y-axis to be aligned with the freestream alpha
.
VortexLattice.Wind
— TypeWind <: AbstractFrame
Reference frame rotated to be aligned with the freestream alpha
and beta
Freestream Parameters
VortexLattice.Freestream
— TypeFreestream([Vinf,] alpha, beta, Omega)
Defines the freestream and rotational velocity properties.
Arguments
Vinf
: Freestream velocityalpha
: angle of attack (rad)beta
: sideslip angle (rad)Omega
: rotation vector (p, q, r) of the body frame about the reference center. Uses standard coordinate system from dynamics (positve p roll right wing down to turn right, positive q is pitch nose up, positive r is yaw nose to the right)
VortexLattice.trajectory_to_freestream
— Functiontrajectory_to_freestream(dt; kwargs...)
Convert trajectory parameters into freestream velocity parameters (see Freestream
) at a collection of time steps.
Arguments:
dt
: Time step vector (seconds)
Keyword Arguments:
Xdot = zeros(length(dt))
: Global frame x-velocity for each time stepYdot = zeros(length(dt))
: Global frame y-velocity for each time stepZdot = zeros(length(dt))
: Global frame z-velocity for each time stepp = zeros(length(dt))
: Angular velocity about x-axis for each time stepq = zeros(length(dt))
: Angular velocity about y-axis for each time stepr = zeros(length(dt))
: Angular velocity about z-axis for each time stepphi0 = 0
: Roll angle for initial time steptheta0 = 0
: Pitch angle for initial time steppsi0 = 0
: Yaw angle for initial time step
Performing an Analysis
VortexLattice.System
— TypeSystem{TF}
Contains pre-allocated storage for internal system variables.
Fields:
AIC
: Aerodynamic influence coefficient matrix from the surface panelsw
: Normal velocity at the control points from external sources and wakesΓ
: Circulation strength of the surface panelsV
: Velocity at the wake vertices for each surfacesurfaces
: Surfaces, represented by matrices of surface panelsproperties
: Surface panel properties for each surfacewakes
: Wake panel properties for each surfacetrefftz
: Trefftz panels associated with each surfacereference
: Pointer to reference parameters associated with the system (seeReference
)freestream
: Pointer to current freestream parameters associated with the system (seeFreestream
)symmetric
: Flags indicating whether each surface is symmetric across the X-Z planenwake
: Number of chordwise wake panels to use from each wake inwakes
,surface_id
: Surface ID for each surface. The finite core model is disabled when calculating the influence of surfaces/wakes that share the same ID.trailing_vortices
: Flags to enable/disable trailing vorticeswake_finite_core
: Flag for each wake indicating whether the finite core model should be enabled when calculating the wake's influence on itself and surfaces/wakes with the same surface ID.near_field_analysis
: Flag indicating whether a near field analysis has been performed for the current system statederivatives
: Flag indicating whether the derivatives with respect to the freestream variables have been calculateddw
: Derivatives of the R.H.S. with respect to the freestream variablesdΓ
: Derivatives of the circulation strength with respect to the freestream variablesdproperties
: Derivatives of the panel properties with respect to the freestream variableswake_shedding_locations
: Wake shedding locations for each surfaceVcp
: Velocity due to surface motion at the control pointsVh
: Velocity due to surface motion at the horizontal bound vortex centersVv
: Velocity due to surface motion at the vertical bound vortex centersVte
: Velocity due to surface motion at the trailing edge verticesdΓdt
: Derivative of the circulation strength with respect to non-dimensional time
VortexLattice.System
— MethodSystem([TF], surfaces; kwargs...)
Return an object of type System
with pre-allocated storage for internal system variables
Arguments:
TF
: Floating point type, defaults to the floating point type used bysurface
surfaces
: Either: - One or more grids of shape (3, nc+1, ns+1) which represents lifting surfaces, or - One or more matrices of surface panels (seeSurfacePanel
) of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panels
Keyword Arguments:
nw
: Number of chordwise wake panels to initialize for each surface. Defaults to zero wake panels for each surface.
VortexLattice.steady_analysis
— Functionsteady_analysis(surfaces, reference, freestream; kwargs...)
Perform a steady vortex lattice method analysis. Return an object of type System
containing the system state.
Arguments
surfaces
:- Vector of grids of shape (3, nc+1, ns+1) which represent lifting surfaces
- Vector of matrices of shape (nc, ns) containing surface panels (see
SurfacePanel
) wherenc
is the number of chordwise panels andns
is the number of spanwise panelsreference
: Reference parameters (seeReference
)freestream
: Freestream parameters (seeFreestream
)
Keyword Arguments
symmetric
: Flag for each surface indicating whether a mirror image across the X-Z plane should be used when calculating induced velocities. Defaults tofalse
for each surfacewakes
: Matrix of wake panels (seeWakePanel
) for each surface. Each matrix has shape (nw, ns) wherenw
is the number of chordwise wake panels andns
is the number of spanwise panels for each surface, defaults to no wake panels for each surfacenwake
: Number of chordwise wake panels to use from each wake inwakes
, defaults to all wake panels for each surfacesurface_id
: Surface ID for each surface. The finite core model is disabled when calculating the influence of surfaces/wakes that share the same ID. By default, all surfaces are assigned their own IDswake_finite_core
: Flag for each wake indicating whether the finite core model should be enabled when calculating a wake's influence on itself and surfaces/wakes with the same surface ID. Defaults totrue
for each surface.trailing_vortices
: Flags to enable/disable trailing vortices for each surface, defaults totrue
for each surfacexhat
: Direction in which to shed trailing vortices, defaults to[1, 0, 0]
additional_velocity
: Function which defines additional velocity as a function of location.fcore
: function which sets the finite core size for each surface based on the chord length and/or the panel width. Defaults to(c, Δs) -> 1e-3
. Only used for grid inputs.calculate_influence_matrix
: Flag indicating whether the aerodynamic influence coefficient matrix has already been calculated. Re-using the same AIC matrix will reduce calculation times when the underlying geometry has not changed. Defaults totrue
. Note that this argument is only valid for the pre-allocated version of this function.near_field_analysis
: Flag indicating whether a near field analysis should be performed to obtain panel velocities, circulation, and forces. Defaults totrue
.derivatives
: Flag indicating whether the derivatives with respect to the freestream variables should be calculated. Defaults totrue
.
VortexLattice.steady_analysis!
— Functionsteady_analysis!(system, surfaces, reference, freestream; kwargs...)
Pre-allocated version of steady_analysis
.
VortexLattice.unsteady_analysis
— Functionunsteady_analysis(surfaces, reference, freestream, dt; kwargs...)
Perform a unsteady vortex lattice method analysis. Return an object of type System
containing the final system state, a matrix of surface panels (see SurfacePanel
for each surface at each time step, a matrix of surface panel properties (see PanelProperties
) for each surface at each time step, and a matrix of wake panels (see WakePanel
) for each surface at each time step.
Arguments
surfaces
:- Grids of shape (3, nc+1, ns+1) which represent lifting surfaces or
- Matrices of surface panels (see
SurfacePanel
) of shape
nc
is the number of chordwise panels andns
is the number of spanwise panels. Alternatively, a vector containing surface shapes/positions at each time step (including att=0
) may be provided to model moving/deforming lifting surfaces.reference
: Reference parameters (seeReference
)freestream
: Freestream parameters for each time step (seeFreestream
)dt
: Time step vector
Keyword Arguments
symmetric
: Flags indicating whether a mirror image (across the X-Z plane) should be used when calculating induced velocities, defaults tofalse
for each surfaceinitial_wakes
: Vector of initial wakes corresponding to each surface, represented by matrices of wake panels (seeWakePanel
) of shape (nw, ns) wherenw
is the number of chordwise wake panels andns
is the number of spanwise panels. Defaults to no wake panels for each surfaceinitial_circulation
: Vector containing the initial circulation of all surface panels in the system. Defaults tozeros(N)
whereN
is the total number of surface panels insurfaces
.nwake
: Maximum number of wake panels in the chordwise direction for each surface. Defaults tolength(dx)
for all surfaces.surface_id
: Surface ID for each surface. The finite core model is disabled when calculating the influence of surfaces/wakes that share the same ID.wake_finite_core
: Flag for each wake indicating whether the finite core model should be enabled when calculating the wake's influence on itself and surfaces/wakes with the same surface ID. Defaults totrue
for each surface.save
: Time indices at which to save the time history, defaults to1:length(dx)
calculate_influence_matrix
: Flag indicating whether the aerodynamic influence coefficient matrix needs to be calculated. Re-using the same AIC matrix will (slightly) reduce calculation times when the underlying geometry has not changed. Defaults totrue
. Note that this argument only affects the pre-allocated version of this function.near_field_analysis
: Flag indicating whether a near field analysis should be performed to obtain panel velocities, circulation, and forces for the final time step. Defaults totrue
.derivatives
: Flag indicating whether the derivatives with respect to the freestream variables should be calculated for the final time step. Defaults totrue
.
VortexLattice.unsteady_analysis!
— Functionunsteady_analysis!(system, surfaces, reference, freestream, dt; kwargs...)
Pre-allocated version of unsteady_analysis
.
Near Field Forces and Moments
VortexLattice.PanelProperties
— TypePanelProperties
Panel specific properties calculated during the vortex lattice method analysis.
Fields
gamma
: Vortex ring circulation strength, normalized by the reference velocityvelocity
: Local velocity at the panel's bound vortex center, normalized by the reference velocitycfb
: Net force on the panel's bound vortex, as calculated using the Kutta-Joukowski theorem, normalized by the reference dynamic pressure and areacfl
: Force on the left bound vortex from this panel's vortex ring, as calculated by the Kutta-Joukowski theorem, normalized by the reference dynamic pressure and areacfr
: Force on the right bound vortex from this panel's vortex ring, as calculated by the Kutta-Joukowski theorem, normalized by the reference dynamic pressure and area
VortexLattice.get_surface_properties
— Functionget_surface_properties(system)
Return a vector of surface panel properties for each surface, stored as matrices of panel properties (see PanelProperties
) of shape (nc, ns) where nc
is the number of chordwise panels and ns
is the number of spanwise panels
VortexLattice.body_forces
— Methodbody_forces(system; kwargs...)
Return the body force coefficients given the panel properties for surfaces
Note that this function assumes that a near-field analysis has already been performed to obtain the panel forces.
Arguments
system
: Object of typeSystem
which holds system properties
Keyword Arguments
frame
: frame in which to returnCF
andCM
, options areBody()
(default),Stability()
, andWind()
`
VortexLattice.body_forces_history
— Functionbody_forces_history(system, surface_history, property_history; frame=Body())
Return the body force coefficients CF
, CM
at each time step in property_history
.
Arguments:
system
: Object of typeSystem
which holds system propertiessurface_history
: Vector of surfaces at each time step, where each surface is represented by a matrix of surface panels (seeSurfacePanel
) of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panelsproperty_history
: Vector of surface properties for each surface at each time step, where surface properties are represented by a matrix of panel properties (seePanelProperties
) of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panels
Keyword Arguments
frame
: frame in which to returnCF
andCM
, options areBody()
(default),Stability()
, andWind()
`
VortexLattice.lifting_line_coefficients
— Functionlifting_line_coefficients(system, r, c; frame=Body())
Return the force and moment coefficients (per unit span) for each spanwise segment of a lifting line representation of the geometry.
This function requires that a near-field analysis has been performed on system
to obtain panel forces.
Arguments
system
: Object of typeSystem
that holds precalculated system properties.r
: Vector with length equal to the number of surfaces, with each element being a matrix with size (3, ns+1) which contains the x, y, and z coordinates of the resulting lifting line coordinatesc
: Vector with length equal to the number of surfaces, with each element being a vector of lengthns+1
which contains the chord lengths at each lifting line coordinate.
Keyword Arguments
frame
: frame in which to returncf
andcm
, possible options areBody()
(default),Stability()
, andWind()
`
Return Arguments:
cf
: Vector with length equal to the number of surfaces, with each element being a matrix with size (3, ns) which contains the x, y, and z direction force coefficients (per unit span) for each spanwise segment.cm
: Vector with length equal to the number of surfaces, with each element being a matrix with size (3, ns) which contains the x, y, and z direction moment coefficients (per unit span) for each spanwise segment.
VortexLattice.lifting_line_coefficients!
— Functionlifting_line_coefficients!(cf, cm, system, r, c; frame=Body())
In-place version of lifting_line_coefficients
Far Field Drag
VortexLattice.far_field_drag
— Functionfar_field_drag(system)
Computes induced drag using the Trefftz plane (far field method).
Note that this function assumes that the circulation distribution has already been computed and is present in system
Arguments
system
: Pre-allocated system properties
far_field_drag(receiving, sending, reference, symmetric)
Computes the induced drag on receiving
from sending
using the Trefftz plane analysis.
Arguments
receiving
: Vector of receiving Trefftz panels (seeTrefftzPanel
)sending
: Vector of sending Trefftz panels (seeTrefftzPanel
)reference
: Reference parameters (seeReference
)symmetric
: Flag indicating whether a mirror image of the panels insurface
, should be used when calculating induced velocities
Body and Stability Derivatives
VortexLattice.body_derivatives
— Functionbody_derivatives(system, surfaces, reference, freestream; kwargs...)
Returns the derivatives of the body forces and moments with respect to the freestream velocity components (u, v, w)
and the angular velocity components (p, q, r)
in the body frame.
The derivatives are returned as two named tuples: dCF, dCM
Note that the derivatives with respect to the freestream variables of the panel forces must have been previously computed and stored in system
Arguments:
system
: Object of typeSystem
which holds system properties
VortexLattice.stability_derivatives
— Functionstability_derivatives(system)
Returns the derivatives of the body forces and moments in the stability frame with respect to the freestream velocity components (alpha, beta)
and the angular velocity components (p, q, r)
in the stability frame.
The derivatives are returned as two named tuples: dCF, dCM
Note that the derivatives with respect to the freestream variables of the panel forces must have been previously computed and stored in system
Arguments:
system
: Object of typeSystem
which holds system properties
Visualization
VortexLattice.write_vtk
— Functionwrite_vtk(name, surfaces, [surface_properties]; kwargs...)
write_vtk(name, wakes; kwargs...)
write_vtk(name, surfaces, wakes, [surface_properties]; kwargs...)
Write geometry from surfaces and/or wakes to Paraview files for visualization.
Arguments
name
: Base name for the generated filessurfaces
:- Vector of grids of shape (3, nc+1, ns+1) which represent lifting surfaces
- Vector of matrices of shape (nc, ns) containing surface panels (see
SurfacePanel
) wherenc
is the number of chordwise panels andns
is the number of spanwise panelswakes
: (optional) Vector of wakes corresponding to each surface, represented by matrices of wake panels (seeWakePanel
) of shape (nw, ns) wherenw
is the number of chordwise wake panels andns
is the number of spanwise panels.surface_properties
: (optional) Vector of surface panel properties for each surface, stored as matrices of panel properties (seePanelProperties
) of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panels
Keyword Arguments:
symmetric
: (required ifsurface_properties
is provided) Flags indicating whether a mirror image (across the X-Z plane) was used when calculating induced velocities for each surface.trailing_vortices
: Flag indicating whether the model uses trailing vortices. Defaults totrue
when wake panels are absent,false
otherwisexhat
: Direction in which trailing vortices extend if used. Defaults to [1, 0, 0].wake_length
: Distance to extend trailing vortices. Defaults to 10metadata
: Dictionary of metadata to include in generated files
write_vtk(name, surface_history, property_history, wake_history; kwargs...)
Writes unsteady simulation geometry to Paraview files for visualization.
Arguments
name
: Base name for the generated filessurface_history
: Vector of surfaces at each time step, where each surface is represented by a matrix of surface panels (seeSurfacePanel
) of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panelsproperty_history
: Vector of surface properties for each surface at each time step, where surface properties are represented by a matrix of panel properties (seePanelProperties
) of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panelswake_history
: Vector of wakes corresponding to each surface at each time step, where each wake is represented by a matrix of wake panels (seeWakePanel
) of shape (nw, ns) wherenw
is the number of chordwise wake panels andns
is the number of spanwise panels.dt
: Time step vector
Keyword Arguments:
symmetric
: (required ifproperties
is provided) Flags indicating whether a mirror image (across the X-Z plane) was used when calculating induced velocities for each surface.wake_length
: Distance to extend trailing vortices. Defaults to 10metadata
: Dictionary of metadata to include in generated files
Private API
Geometry
VortexLattice.linearinterp
— Functionlinearinterp(eta, rstart, rend)
Linearly interpolate between rstart
and rend
where eta is the fraction between 0 (rstart) and 1 (rend)
VortexLattice.spanwise_spacing
— Functionspanwise_spacing(n, spacing::AbstractSpacing)
Distribute n
panel endpoints and n-1
panel midpoints on the interval between 0 and 1 according to the discretization strategy in spacing
.
VortexLattice.chordwise_spacing
— Functionchordwise_spacing(n, spacing::AbstractSpacing)
Distribute n
panel edge, n-1
vortex, and n-1
control point chordwise locations on the interval between 0 and 1 according to the discretization strategy in spacing
.
VortexLattice.interpolate_grid
— Functioninterpolate_grid(xyz, eta, interp; xdir=0, ydir=1)
Interpolates the grid xyz
along direction dir
Arguments
xyz
: Grid of size (3, ni, nj)eta
: New (normalized) coordinates in directiondir
(0 <= eta <= 1)interp
: Interpolation method of formypt = f(x,y,xpt)
xdir
: Independent variable direction, defaults to arc lengthydir
: Dependent variable directionxyz
(i=1
,j=2
)
VortexLattice.update_surface_panels!
— Functionupdate_surface_panels!(surface, grid; fcore = (c, Δs) -> 1e-3)
Updates the surface panels in surface
to correspond to the grid coordinates in grid
.
VortexLattice.trailing_edge_points
— Functiontrailing_edge_points(surface[s])
Return points on the trailing edge of each surface.
VortexLattice.repeated_trailing_edge_points
— Functionrepeated_trailing_edge_points(surface[s])
Generates a dictionary of the form Dict((isurf, i) => [(jsurf1, j1), (jsurf2, j2)...]
which defines repeated trailing edge points. Trailing edge point i
on surface isurf
is repeated on surface jsurf1
at point j1
, jsurf2
at point j2
, and so forth.
VortexLattice.flipy
— Functionflipy(r)
Flip sign of y-component of vector r
(used for symmetry)
VortexLattice.on_symmetry_plane
— Functionon_symmetry_plane(args...; tol=eps())
Test whether points in args
are on the symmetry plane (y = 0)
VortexLattice.not_on_symmetry_plane
— Functionnot_on_symmetry_plane(args...; tol=eps())
Test whether and of the points in args
are not on the symmetry plane (y = 0)
Surface Panels
VortexLattice.top_left
— Functiontop_left(panel::SurfacePanel)
Return the top left vertex of the vortex ring associated with panel
VortexLattice.top_center
— Functiontop_center(panel::SurfacePanel)
Return the top center vertex of the vortex ring associated with panel
VortexLattice.top_right
— Functiontop_right(panel::SurfacePanel)
Return the top right vertex of the vortex ring associated with panel
VortexLattice.bottom_left
— Functionbottom_left(panel::SurfacePanel)
Return the bottom left vertex of the vortex ring associated with panel
VortexLattice.bottom_center
— Functionbottom_center(panel::SurfacePanel)
Return the bottom center vertex of the vortex ring associated with panel
VortexLattice.bottom_right
— Functionbottom_right(panel::SurfacePanel)
Return the bottom right vertex of the vortex ring associated with panel
VortexLattice.controlpoint
— Functioncontrolpoint(panel::SurfacePanel)
Return the control point of panel
(typically located at the 3/4 chord)
VortexLattice.normal
— Methodnormal(panel::SurfacePanel)
Return the normal vector of panel
at the panel control point
VortexLattice.get_core_size
— Functionget_core_size(panel::SurfacePanel)
Return the core size (smoothing parameter) corresponding to the vortex ring associated with panel
VortexLattice.reflect
— Methodreflect(panel::SurfacePanel)
Reflect panel
across the X-Z plane.
VortexLattice.left_center
— Functionleft_center(panel::SurfacePanel)
Return the center of the left bound vortex of the vortex ring associated with panel
VortexLattice.right_center
— Functionright_center(panel::SurfacePanel)
Return the center of the right bound vortex of the vortex ring associated with panel
VortexLattice.top_vector
— Functiontop_vector(panel::SurfacePanel)
Return the path of the top bound vortex of the vortex ring associated with panel
VortexLattice.left_vector
— Functionleft_vector(panel)
Return the path of the left bound vortex of the vortex ring associated with panel
VortexLattice.right_vector
— Functionright_vector(panel)
Return the path of the right bound vortex of the vortex ring associated with panel
VortexLattice.bottom_vector
— Functionbottom_vector(panel)
Return the path of the bottom bound vortex of the vortex ring associated with panel
Wake Panels
VortexLattice.update_wake_shedding_locations!
— Functionupdate_wake_shedding_locations!(wakes, wake_shedding_locations,
surfaces, ref, fs, dt, additional_velocity, Vte, nwake, eta)
Update the wake shedding locations. Also update the first chordwise wake panels to account for the new wake shedding location
Arguments
wakes
: Vector of wakes corresponding to each surface, represented by matrices of wake panels (seeWakePanel
) of shape (nw, ns) wherenw
is the number of chordwise wake panels andns
is the number of spanwise panels.wake_shedding_locations
: Shedding location coordinates for each surface for each trailing edge vertex.surfaces
: Vector of surfaces, represented by matrices of surface panels (seeSurfacePanel
of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panelsreference
: Reference parameters (seeReference
)freestream
: Freestream parameters (seeFreestream
)dt
: Time step (seconds)additional_velocity
: Function defining additional velocity fieldVte
: Velocity experienced at the trailing edge due to surface motion.nwake
: Number of chordwise wake panels to use from each wake inwakes
eta
: Time step fraction used to define separation between trailing edge and wake shedding location. Typical values range from 0.2-0.3.
VortexLattice.set_circulation_strength
— Functionset_circulation_strength(panel::WakePanel, gamma)
Return a copy of panel
with the circulation strength gamma
VortexLattice.circulation_strength
— Functioncirculation_strength(panel::WakePanel)
Return the circulation strength of the wake panel.
VortexLattice.get_wake_velocities!
— Functionget_wake_velocities!(wake_velocities, surfaces, wakes, ref, fs, Γ,
additional_velocity, Vte, symmetric, repeated_points, nwake,
surface_id, wake_finite_core, wake_shedding_locations, trailing_vortices, xhat)
Arguments
wake_velocities
: Velocities at the corners of the wake panels inwakes
surfaces
: Vector of surfaces, represented by matrices of surface panels (seeSurfacePanel
of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panelswakes
: Vector of wakes corresponding to each surface, represented by matrices of wake panels (seeWakePanel
) of shape (nw, ns) wherenw
is the number of chordwise wake panels andns
is the number of spanwise panels.reference
: Reference parameters (seeReference
)freestream
: Freestream parameters (seeFreestream
)Γ
: Circulation of all surface panels stored in a single vectoradditional_velocity
: Function defining additional velocity fieldVte
: Velocity at the trailing edge vertices on each surface due to surface motionsymmetric
: (required) Flag for each surface indicating whether a mirror image across the X-Z plane should be used when calculating induced velocitiesrepeated_points
: Dictionary of the formDict((isurf, i) => [(jsurf1, j1), (jsurf2, j2)...]
which defines repeated trailing edge points. Trailing edge pointi
on surfaceisurf
is repeated on surfacejsurf1
at pointj1
,jsurf2
at pointj2
, and so forth. Seerepeated_trailing_edge_points
nwake
: Number of chordwise wake panels to use from each wake inwakes
, defaults to all provided wake panelssurface_id
: Surface ID for each surface. The finite core model is disabled when calculating the influence of surfaces/wakes that share the same ID.wake_finite_core
: Flag for each wake indicating whether the finite core model should be enabled when calculating the wake's influence on itself and surfaces/wakes with the same surface ID. Defaults totrue
for each surface.wake_shedding_locations
: Shedding location coordinates for each surface for each trailing edge vertex.trailing_vortices
: Flags to enable/disable trailing vortices, defaults totrue
for each surfacexhat
: Direction in which to shed trailing vortices, defaults to [1, 0, 0]
VortexLattice.translate_wake
— Functiontranslate_wake(panel, wake_velocities, dt)
Return a translated copy of the wake panel panel
given the wake corner velocities wake_velocities
and the time step dt
Arguments
panel
: Wake panel (of typeWakePanel
)wake_velocities
: Matrix containing the velocities at each of the four corners ofpanel
dt
: Time step (seconds)
VortexLattice.translate_wake!
— Functiontranslate_wake!(wake, wake_velocities, dt; nwake = size(wake, 1))
Translate the wake panels in wake
given the corner velocities wake_velocities
and the time step dt
.
Arguments
wake
: Matrix of wake panels (seeWakePanel
) of shape (nw, ns) wherenw
is the number of chordwise wake panels andns
is the number of spanwise panels, defaults to no wake panelswake_velocities
: Velocities at each of the vertices corresponding to the wake panels inwake
dt
: Time step
Keyword Arguments
nwake
: Number of chordwise wake panels to use fromwake
, defaults to all provided wake panels
VortexLattice.shed_wake!
— Functionshed_wake!(wake, wake_shedding_locations, wake_velocities, dt, surface, Γ, nwake)
Shed a new wake panel from the wake shedding locations and translate existing wake panels.
Arguments
wake
: Matrix of wake panels (seeWakePanel
) of shape (nw, ns) wherenw
is the number of chordwise wake panels andns
is the number of spanwise panelswake_shedding_locations
: Vector of lengthns
which stores the coordinates where wake panels are shed from the trailing edge ofsurface
.wake_velocities
: Velocities at each of the vertices corresponding to the wake panels inwake
dt
: Time step (seconds)surface
: Matrix of surface panels (seeSurfacePanel
) of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panelsΓ
: Circulation strength of each surface panel insurface
nwake
: Number of chordwise wake panels to use fromwake
, defaults to all provided wake panels
shed_wake!(wakes, wake_shedding_locations, wake_velocities, dt, surfaces, Γ, nwake)
Shed a new wake panel from the wake shedding locations and translate existing wake panels.
Arguments
wakes
: Vector of wakes corresponding to each surface, represented by matrices of wake panels (seeWakePanel
) of shape (nw, ns) wherenw
is the number of chordwise wake panels andns
is the number of spanwise panels.wake_shedding_locations
: Shedding location coordinates for each surface for each trailing edge vertex.wake_velocities
: Velocities at each of the vertices corresponding to the wake panels inwake
dt
: Time step (seconds)surfaces
: Vector of surfaces, represented by matrices of surface panels (seeSurfacePanel
of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panelsΓ
: Circulation strength of each surface panel insurfaces
nwake
: Number of chordwise wake panels to use from each wake inwakes
, defaults to all provided wake panels
VortexLattice.rowshift!
— Functionrowshift!(A)
Circularly shifts the rows of a matrix down one row.
Induced Velocity
VortexLattice.bound_induced_velocity
— Functionbound_induced_velocity(r1, r2, finite_core, core_size)
Compute the induced velocity (per unit circulation) for a bound vortex, at a control point located at r1
relative to the start of the bound vortex and r2
relative to the end of the bound vortex
VortexLattice.trailing_induced_velocity
— Functiontrailing_induced_velocity(r1, r2, xhat, finite_core, core_size)
Compute the induced velocity (per unit circulation) for a vortex trailing in the xhat
direction, at a control point located at r
relative to the start of the trailing vortex.
VortexLattice.ring_induced_velocity
— Functionring_induced_velocity(rcp, r11, r12, r21, r22; finite_core = false,
core_size = 0.0, symmetric = false, xhat = [1,0,0], top = true, bottom = true,
left = true, right = true, left_trailing = false, right_trailing = false,
reflected_top = true, reflected_bottom = true, reflected_left = true,
reflected_right = true, reflected_left_trailing = false, reflected_right_trailing = false)
Compute the induced velocity (per unit circulation) for a vortex ring defined by the corners r11
, r12
, r21
, and r22
at a control point located at rcp
Also returns the induced velocity resulting from shared edges with panels on the top, bottom, left, and right sides of the panel described by r11
, r12
, r21
, and r22
.
ring_induced_velocity(rcp, panel; kwargs...)
Compute the velocity (per unit circulation) induced by panel
at a control point located at rcp
Also returns the velocity induced by the shared edges of adjacent panels on the top, bottom, left, and right sides of panel
.
VortexLattice.influence_coefficients!
— Functioninfluence_coefficients!(AIC, surface; kwargs...)
Construct the aerodynamic influence coefficient matrix for a single surface.
Arguments:
surface
: Matrix of surface panels (seeSurfacePanel
) of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panels
Keyword Arguments
symmetric
: Flag indicating whether a mirror image of the panels insurface
should be used when calculating induced velocities.wake_shedding_locations
: Wake shedding locations for the trailing edge panels insurface
trailing_vortices
: Flag to enable/disable trailing vortices. Defaults totrue
.xhat
: Direction in which trailing vortices are shed iftrailing_vortices = true
. Defaults to [1, 0, 0]
influence_coefficients!(AIC, surfaces; kwargs...)
Construct the aerodynamic influence coefficient matrix for multiple surfaces.
Arguments:
surfaces
: Vector of surfaces, represented by matrices of surface panels (seeSurfacePanel
of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panels
Keyword Arguments:
symmetric
: Flags indicating whether a mirror image (across the X-Z plane) should be used when calculating induced velocities. Defaults tofalse
for each surface.surface_id
: ID for each surface. May be used to deactivate the finite core model by setting all surface ID's to the same value. Defaults to a unique ID for each surfacewake_shedding_locations
: Wake shedding locations for the trailing edge panels of each surface insurfaces
trailing_vortices
: Flags to indicate whether trailing vortices are used for each surface. Defaults totrue
for each surface.xhat
: Direction in which trailing vortices are shed iftrailing_vortices = true
. Defaults to [1, 0, 0]
influence_coefficients!(AIC, receiving, sending; kwargs...)
Compute the AIC coefficients corresponding to the influence of the panels in sending
on the panels in receiving
.
Keyword Arguments
finite_core
: Flag indicating whether the finite core model is enabled. Defaults totrue
symmetric
: Flag indicating whether sending panels should be mirrored across the X-Z planewake_shedding_locations
: Wake shedding locations for the trailing edge panels insending
trailing_vortices
: Indicates whether trailing vortices are used. Defaults totrue
.xhat
: Direction in which trailing vortices are shed iftrailing_vortices = true
. Defaults to [1, 0, 0]
VortexLattice.update_trailing_edge_coefficients!
— Functionupdate_trailing_edge_coefficients!(AIC, surface; kwargs...)
Construct the aerodynamic influence coefficient matrix for a single surface.
Arguments:
surface
: Matrix of panels of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panels
Keyword Arguments
symmetric
: Flag indicating whether a mirror image of the panels insurface
should be used when calculating induced velocities.trailing_vortices
: Flag to enable/disable trailing vortices.xhat
: Direction in which trailing vortices are shed iftrailing_vortices = true
. Defaults to [1, 0, 0]
update_trailing_edge_coefficients!(AIC, surfaces; kwargs...)
Construct the aerodynamic influence coefficient matrix for multiple surfaces.
Arguments:
surfaces
: Vector of surfaces, represented by matrices of panels of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panels
Keyword Arguments:
symmetric
: Flags indicating whether a mirror image (across the X-Z plane) should be used when calculating induced velocities. Defaults tofalse
for each surface.wake_shedding_locations
: Shedding location coordinates for each surface for each trailing edge vertex.surface_id
: ID for each surface. May be used to deactivate the finite core model by setting all surface ID's to the same value. Defaults to a unique ID for each surfacetrailing_vortices
: Flags to indicate whether trailing vortices are used for each surface. Defaults totrue
for each surface.xhat
: Direction in which trailing vortices are shed iftrailing_vortices = true
. Defaults to [1, 0, 0]
update_trailing_edge_coefficients!(AIC, receiving, sending; kwargs...)
Update the AIC coefficients corresponding to the influence of the trailing edge panels in sending
on the panels in receiving
.
Keyword Arguments
finite_core
: Flag indicating whether the finite core model is enabled. Defaults totrue
symmetric
: Flag indicating whether sending panels should be mirrored across the X-Z planewake_shedding_locations
: Wake shedding locations for the trailing edge panels insending
trailing_vortices
: Indicates whether trailing vortices are used. Defaults totrue
.xhat
: Direction in which trailing vortices are shed iftrailing_vortices = true
. Defaults to [1, 0, 0]
VortexLattice.induced_velocity
— Functioninduced_velocity(rcp, surface, Γ; kwargs...)
Compute the velocity induced by the grid of panels in surface
at control point rcp
induced_velocity(I::CartesianIndex, surface, Γ; kwargs...)
Compute the velocity induced by the grid of panels in surface
on the top bound vortex of panel I
in surface
.
induced_velocity(is::Integer, surface, Γ; kwargs...)
Compute the velocity induced by the grid of panels in surface
at the trailing edge vertex corresponding to index is
induced_velocity(rcp, wake::AbstractMatrix{<:WakePanel}; kwargs...)
Compute the velocity induced by the grid of wake panels in wake
at control point rcp
induced_velocity(I::CartesianIndex, wake; kwargs...)
Compute the induced velocity from the grid of wake panels in wake
at the vertex corresponding to index I
induced_velocity(rcp, surface, Γ = nothing, dΓ = nothing; kwargs...)
Compute the induced velocity from the grid of panels in surface
at rcp
using the circulation strengths provided in Γ.
Keyword Arguments
nc
: Number of panels in the chordwise direction. Defaults tosize(surface, 1)
ns
: Number of panels in the spanwise direction. Defaults tosize(surface, 2)
finite_core
: Flag indicating whether the finite core model should be usedsymmetric
: Flag indicating whether sending panels should be mirrored across the X-Z planewake_shedding_locations
: Wake shedding locations for the trailing edge panels insurface
trailing_vortices
: Indicates whether trailing vortices are used. Defaults totrue
.xhat
: Direction in which trailing vortices are shed iftrailing_vortices = true
. Defaults to [1, 0, 0]skip_leading_edge = false
: Indicates whether to skip the leading edge. This flag may be used to skip calculating the leading bound vortex of a wake when its influence cancels exactly with the trailing bound vortex of a surface.skip_inside_edges = false
: Indicates whether to skip all horizontal bound vortices except those located at the leading and trailing edges. This flag may be used to skip calculating a wake's (internal) horizontal bound vortices during steady state simulations since the influence of adjacent wake panels in a chordwise strip cancels exactly in steady state simulations.skip_trailing_edge = false
: Indicates whether to skip the trailing edge. The trailing edge is always skipped iftrailing_vortices = true
skip_top
: Tuple containing panel indices whose top bound vortex is coincident withrcp
and should therefore be skipped.skip_bottom
: Tuple containing panel indices whose bottom bound vortex is coincident withrcp
and should therefore be skipped.skip_left
: Tuple containing panel indices whose left bound vortex is coincident withrcp
and should therefore be skipped.skip_right
: Tuple containing panel indices whose right bound vortex is coincident withrcp
and should therefore be skipped.skip_left_trailing
: Tuple containing panel indices whose left trailing vortex is coincident withrcp
and should therefore be skipped.skip_right_trailing
: Tuple containing panel indices whose right trailing vortex is coincident withrcp
and should therefore be skipped.
VortexLattice.induced_velocity_derivatives
— Functioninduced_velocity_derivatives(rcp, surface, Γ, dΓ; kwargs...)
Compute the velocity induced by the grid of panels in surface
at control point rcp
and its derivatives with respect to the freestream variables
induced_velocity_derivatives(I::CartesianIndex, surface, Γ, dΓ; kwargs...)
Compute the velocity induced by the grid of panels in surface
on the top bound of panel I
in surface
and its derivatives with respect to the freestream variables
Freestream
VortexLattice.body_to_stability
— Functionbody_to_stability(freestream)
Construct a rotation matrix from the body axis to the stability axis.
VortexLattice.body_to_wind
— Functionbody_to_wind(freestream)
Construct a rotation matrix from the body axis to the wind axis
VortexLattice.stability_to_body
— Functionstability_to_body(freestream)
Construct a rotation matrix from the stability axis to the body axis
VortexLattice.stability_to_wind
— Functionstability_to_wind(freestream)
Construct a rotation matrix from the stability axis to the wind axis
VortexLattice.wind_to_body
— Functionwind_to_body(freestream)
Construct a rotation matrix from the wind axis to the body axis
VortexLattice.wind_to_stability
— Functionwind_to_stability(freestream)
Construct a rotation matrix from the wind axis to the stability axis
VortexLattice.body_to_stability_alpha
— Functionbody_to_stability_alpha(freestream)
Construct a rotation matrix from the body axis to the stability axis and its derivative with respect to alpha
VortexLattice.body_to_wind_derivatives
— Functionbody_to_wind_derivatives(freestream)
Construct a rotation matrix from the body axis to the wind axis and its derivatives with respect to alpha
and beta
VortexLattice.stability_to_body_alpha
— Functionstability_to_body(freestream)
Construct a rotation matrix from the stability axis to the body axis and its derivative with respect to alpha
VortexLattice.stability_to_wind_beta
— Functionstability_to_wind_beta(freestream)
Construct a rotation matrix from the stability axis to the wind axis and its derivative with respect to beta
VortexLattice.wind_to_body_derivatives
— Functionwind_to_body_derivatives(freestream)
Construct a rotation matrix from the wind axis to the body axis and its derivatives with respect to alpha
and beta
VortexLattice.wind_to_stability_beta
— Functionwind_to_stability_beta(freestream)
Construct a rotation matrix from the wind axis to the stability axis and its derivative with respect to beta
VortexLattice.freestream_velocity
— Functionfreestream_velocity(freestream)
Computes the freestream velocity
VortexLattice.freestream_velocity_derivatives
— Functionfreestream_velocity_derivatives(freestream)
Computes the freestream velocity
VortexLattice.rotational_velocity
— Functionrotational_velocity(r, freestream, reference)
Compute the velocity due to body rotations about the reference center
VortexLattice.rotational_velocity_derivatives
— Functionrotational_velocity_derivatives(r, freestream, reference)
Compute the velocity due to body rotations about the reference center and its derivatives with respect to (p, q, r)
VortexLattice.get_surface_velocities!
— Functionget_surface_velocities!(Vcp, Vh, Vv, Vte, current_surface, previous_surface, dt)
Calculate the velocities experienced by the surface at the control points, horizontal bound vortex centers, vertical bound vortex centers and trailing edge vertices due to surface motion.
Circulation
VortexLattice.normal_velocity!
— Functionnormal_velocity!(w, surfaces, wakes, ref, fs; additional_velocity,
Vcp, symmetric, nwake, surface_id, wake_finite_core, trailing_vortices, xhat)
Compute the downwash at the control points on surfaces
due to the freestream velocity, rotational velocity, additional velocity field, surface motion, and induced velocity from the wake panels.
This forms the right hand side of the circulation linear system solve.
VortexLattice.normal_velocity_derivatives!
— Functionnormal_velocity_derivatives!(w, dw, surfaces, wakes, ref, fs;
additional_velocity, Vcp, symmetric, nwake, surface_id,
wake_finite_core, trailing_vortices, xhat)
Compute the downwash at the control points on surfaces
due to the freestream velocity, rotational velocity, additional velocity field, surface motion, and induced velocity from the wake panels. Also calculate its derivatives with respect to the freestream parameters.
This forms the right hand side of the circulation linear system solve (and its derivatives).
VortexLattice.circulation
— Functioncirculation(AIC, w)
Solve for the circulation distribution.
VortexLattice.circulation!
— Functioncirculation!(Γ, AIC, w)
Pre-allocated version of circulation
VortexLattice.circulation_derivatives
— Functioncirculation_derivatives(AIC, w, dw)
Solve for the circulation distribution and its derivatives with respect to the freestream parameters.
VortexLattice.circulation_derivatives!
— Functioncirculation_derivatives!(Γ, dΓ, AIC, w, dw)
Pre-allocated version of circulation_derivatives
Time-Domain Analysis
VortexLattice.propagate_system!
— Functionpropagate_system!(system, [surfaces, ] freestream, dt; kwargs...)
Propagate the state variables in system
forward one time step using the unsteady vortex lattice method system of equations.
Arguments
system
: Object of typesystem
which contains the current system statesurfaces
: Surface locations at the end of this time step. If omitted, surfaces are assumed to be stationary.freestream
: Freestream parameters corresponding to this time step.dt
: Time increment
Keyword Arguments
additional_velocity
: Function which defines additional velocity as a function of location.repeated_points
: Dictionary of the formDict((isurf, i) => [(jsurf1, j1), (jsurf2, j2)...]
which defines repeated trailing edge points. Trailing edge pointi
on surfaceisurf
is repeated on surfacejsurf1
at pointj1
,jsurf2
at pointj2
, and so forth. Seerepeated_trailing_edge_points
nwake
: Number of wake panels in the chordwise direction for each surface.eta
: Time step fraction used to define separation between trailing edge and wake shedding location. Typical values range from 0.2-0.3.calculate_influence_matrix
: Flag indicating whether the aerodynamic influence coefficient matrix needs to be calculated. If argumentsurfaces
is provided the influence matrix will always be recalculated.near_field_analysis
: Flag indicating whether a near field analysis should be performed to obtain panel velocities, circulation, and forces.derivatives
: Flag indicating whether the derivatives with respect to the freestream variables should be calculated.
Near-Field Analysis
VortexLattice.near_field_forces!
— Functionnear_field_forces!(properties, surfaces, wakes, reference, freestream, Γ;
dΓdt, additional_velocity, Vh, Vv, symmetric, nwake, surface_id,
wake_finite_core, wake_shedding_locations, trailing_vortices, xhat)
Calculate local panel forces in the body frame.
VortexLattice.near_field_forces_derivatives!
— Functionnear_field_forces_derivatives!(properties, dproperties, surfaces, reference,
freestream, Γ, dΓ; dΓdt, additional_velocity, Vh, Vv, symmetric, nwake,
surface_id, wake_finite_core, wake_shedding_locations, trailing_vortices, xhat)
Version of near_field_forces!
that also calculates the derivatives of the local panel forces with respect to the freestream variables.
VortexLattice.body_forces
— Methodbody_forces(surfaces, properties, reference, freestream, symmetric; kwargs...)
Return the body force coefficients given the panel properties for surfaces
Note that this function assumes that a near-field analysis has already been performed to obtain the panel forces.
Arguments:
surfaces
: Collection of surfaces, where each surface is represented by a matrix of surface panels (seeSurfacePanel
) of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panelsproperties
: Surface properties for each surface, where surface properties for each surface are represented by a matrix of panel properties (seePanelProperties
) of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panelsreference
: Reference parameters (seeReference
)freestream
: Freestream parameters (see [Freestream
]@ref)symmetric
: (required) Flag for each surface indicating whether a mirror image (across the X-Z plane) was used when calculating induced velocitiesframe
: frame in which to returnCF
andCM
, options areBody()
(default),Stability()
, andWind()
VortexLattice.body_forces_derivatives
— Functionbody_forces_derivatives(system)
Return the body force coefficients for the system
and their derivatives with respect to the freestream variables
Note that this function assumes that a near-field analysis has already been performed to obtain the panel forces.
Arguments:
system
: Object of typeSystem
which holds system properties
VortexLattice.body_to_frame
— Functionbody_to_frame(CF, CM, reference, freestream, frame)
Transform the coefficients CF
and CM
from the body frame to the frame specified in frame
Far-Field
VortexLattice.TrefftzPanel
— TypeTrefftzPanel{TF}
Panel in the Trefftz plane.
VortexLattice.normal
— Methodnormal(panel::TrefftzPanel)
Return the normal vector of panel
, including magnitude
VortexLattice.trefftz_panels
— Functiontrefftz_panels(surface[s], freestream, Γ)
Constructs a set of panels for Trefftz plane calculations
VortexLattice.trefftz_panel_induced_drag
— Functiontrefftz_panel_induced_drag(receiving::TrefftzPanel, sending::TrefftzPanel; kwargs...)
Induced drag on receiving
panel induced by sending
panel.
Keyword Arguments
symmetric
: Flag indicating whether a mirror image ofsending
should be used when calculating the induced drag
VortexLattice.vortex_induced_drag
— Functionvortex_induced_drag(rj, Γj, ri, Γi, ni)
Return induced drag from vortex j
induced on panel i
Visualization
VortexLattice.write_vtk!
— Functionwrite_vtk!(vtmfile, surface, [surface_properties]; kwargs...)
Writes geometry to Paraview files for visualization.
Arguments
vtmfile
: Multiblock file handlesurface
: Matrix of surface panels (seeSurfacePanel
) of shape (nc, ns) wherenc
is the number of chordwise panels andns
is the number of spanwise panelssurface_properties
: (optional) Matrix of panel properties for each non-wake panel where each element of the matrix is of typePanelProperties
.
Keyword Arguments:
symmetric
: (required ifproperties
is provided) Flag indicating whether a mirror image (across the X-Z plane) was used when calculating induced velocities.trailing_vortices = true
: Flag indicating whether the model uses trailing vorticesxhat = [1, 0, 0]
: Direction in which trailing vortices extend if usedwake_length = 10
: Distance to extend trailing vorticeswake_circulation = zeros(size(surfaces, 2))
: Contribution to the trailing edge circulation from the wake attached to this surfacemetadata = Dict()
: Dictionary of metadata to include in generated files
write_vtk!(vtmfile, wake; kwargs...)
Writes geometry to Paraview files for visualization.
Arguments
vtmfile
: Paraview file handlewake
: Matrix of wake panels (seeWakePanel
) of shape (nw, ns) wherenw
is the number of chordwise wake panels andns
is the number of spanwise panels
Keyword Arguments:
symmetric
: (required) Flag indicating whether a mirror image (across the X-Z plane) was used when calculating induced velocities.trailing_vortices = false
: Flag indicating whether the model uses trailing vorticesxhat = [1, 0, 0]
: Direction in which trailing vortices extend if usedwake_length = 10
: Distance to extend trailing vorticessurface_circulation = zeros(size(wake, 2))
: Contribution to the leading edge circulation from the surface attached to this wake.metadata = Dict()
: Dictionary of metadata to include in generated files
Index
VortexLattice.AbstractFrame
VortexLattice.AbstractSpacing
VortexLattice.Body
VortexLattice.Cosine
VortexLattice.Freestream
VortexLattice.PanelProperties
VortexLattice.Reference
VortexLattice.Sine
VortexLattice.Stability
VortexLattice.SurfacePanel
VortexLattice.SurfacePanel
VortexLattice.System
VortexLattice.System
VortexLattice.TrefftzPanel
VortexLattice.Uniform
VortexLattice.WakePanel
VortexLattice.WakePanel
VortexLattice.Wind
VortexLattice.body_derivatives
VortexLattice.body_forces
VortexLattice.body_forces
VortexLattice.body_forces_derivatives
VortexLattice.body_forces_history
VortexLattice.body_to_frame
VortexLattice.body_to_stability
VortexLattice.body_to_stability_alpha
VortexLattice.body_to_wind
VortexLattice.body_to_wind_derivatives
VortexLattice.bottom_center
VortexLattice.bottom_left
VortexLattice.bottom_right
VortexLattice.bottom_vector
VortexLattice.bound_induced_velocity
VortexLattice.chordwise_spacing
VortexLattice.circulation
VortexLattice.circulation!
VortexLattice.circulation_derivatives
VortexLattice.circulation_derivatives!
VortexLattice.circulation_strength
VortexLattice.controlpoint
VortexLattice.far_field_drag
VortexLattice.flipy
VortexLattice.freestream_velocity
VortexLattice.freestream_velocity_derivatives
VortexLattice.get_core_size
VortexLattice.get_surface_properties
VortexLattice.get_surface_velocities!
VortexLattice.get_wake_velocities!
VortexLattice.grid_to_surface_panels
VortexLattice.import_vsp
VortexLattice.induced_velocity
VortexLattice.induced_velocity_derivatives
VortexLattice.influence_coefficients!
VortexLattice.interpolate_grid
VortexLattice.left_center
VortexLattice.left_vector
VortexLattice.lifting_line_coefficients
VortexLattice.lifting_line_coefficients!
VortexLattice.lifting_line_geometry
VortexLattice.lifting_line_geometry!
VortexLattice.linearinterp
VortexLattice.near_field_forces!
VortexLattice.near_field_forces_derivatives!
VortexLattice.normal
VortexLattice.normal
VortexLattice.normal_velocity!
VortexLattice.normal_velocity_derivatives!
VortexLattice.not_on_symmetry_plane
VortexLattice.on_symmetry_plane
VortexLattice.propagate_system!
VortexLattice.read_degengeom
VortexLattice.reflect
VortexLattice.reflect
VortexLattice.repeated_trailing_edge_points
VortexLattice.right_center
VortexLattice.right_vector
VortexLattice.ring_induced_velocity
VortexLattice.rotate
VortexLattice.rotate!
VortexLattice.rotational_velocity
VortexLattice.rotational_velocity_derivatives
VortexLattice.rowshift!
VortexLattice.set_circulation_strength
VortexLattice.set_normal
VortexLattice.shed_wake!
VortexLattice.spanwise_spacing
VortexLattice.stability_derivatives
VortexLattice.stability_to_body
VortexLattice.stability_to_body_alpha
VortexLattice.stability_to_wind
VortexLattice.stability_to_wind_beta
VortexLattice.steady_analysis
VortexLattice.steady_analysis!
VortexLattice.top_center
VortexLattice.top_left
VortexLattice.top_right
VortexLattice.top_vector
VortexLattice.trailing_edge_points
VortexLattice.trailing_induced_velocity
VortexLattice.trajectory_to_freestream
VortexLattice.translate
VortexLattice.translate!
VortexLattice.translate_wake
VortexLattice.translate_wake!
VortexLattice.trefftz_panel_induced_drag
VortexLattice.trefftz_panels
VortexLattice.unsteady_analysis
VortexLattice.unsteady_analysis!
VortexLattice.update_surface_panels!
VortexLattice.update_trailing_edge_coefficients!
VortexLattice.update_wake_shedding_locations!
VortexLattice.vortex_induced_drag
VortexLattice.wind_to_body
VortexLattice.wind_to_body_derivatives
VortexLattice.wind_to_stability
VortexLattice.wind_to_stability_beta
VortexLattice.wing_to_surface_panels
VortexLattice.write_vtk
VortexLattice.write_vtk!