Library

Public API

Generating Lifting Surfaces

VortexLattice.SineType
Sine()

Sine-spaced discretization scheme. Using sine-spacing on the right half of a wing effectively results in cosine spacing once symmetry is applied.

source
VortexLattice.CosineType
Cosine()

Cosine-spaced discretization scheme. This is typically one of the most accurate spacing schemes for spanwise spacing.

source
VortexLattice.SurfacePanelType
SurfacePanel{TF}

Lifting surface panel with attached vortex ring

Fields

  • rtl: position of the left side of the top bound vortex
  • rtc: position of the center of the top bound vortex
  • rtr: position of the right side of the top bound vortex
  • rbl: position of the left side of the bottom bound vortex
  • rbc: position of the center of the bottom bound vortex
  • rbr: position of the right side of the bottom bound vortex
  • rcp: position of the panel control point
  • ncp: normal vector at the panel control point
  • core_size: finite core size (for use when the finite core smoothing model is enabled)
  • chord: panel chord length (for determining unsteady forces)
source
VortexLattice.SurfacePanelMethod
SurfacePanel(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 vortex
  • rtr: position of the right side of the top bound vortex
  • rbl: position of the left side of the bottom bound vortex
  • rbr: position of the right side of the bottom bound vortex
  • rcp: position of the panel control point
  • ncp: normal vector at the panel control point
  • core_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
source
VortexLattice.WakePanelType
WakePanel{TF}

SurfacePanel used for modeling wakes.

Fields

  • rtl: position of the left side of the top bound vortex
  • rtr: position of the right side of the top bound vortex
  • rbl: position of the left side of the bottom bound vortex
  • rbr: position of the right side of the bottom bound vortex
  • core_size: finite core size (for use when the finite core smoothing model is enabled)
  • gamma: circulation strength of the panel
source
VortexLattice.WakePanelMethod
WakePanel(rtl, rtr, rbl, rbr, core_size, gamma)

Construct and return a wake panel.

Arguments

  • rtl: position of the left side of the top bound vortex
  • rtr: position of the right side of the top bound vortex
  • rbl: position of the left side of the bottom bound vortex
  • rbr: position of the right side of the bottom bound vortex
  • core_size: finite core size
  • gamma: circulation strength of the panel
source
VortexLattice.grid_to_surface_panelsFunction
grid_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 to false.
  • 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
source
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 where i corresponds to the chordwise direction and j corresponds to the spanwise direction.
  • ns: number of spanwise panels
  • nc: number of chordwise panels

Keyword Arguments

  • mirror: mirror the geometry across the X-Z plane? defaults to false.
  • 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 to Cosine()
  • spacing_c: chordwise discretization scheme, defaults to Uniform()
  • interp_s: spanwise interpolation function, defaults to linear interpolation
  • interp_c: chordwise interpolation function, defaults to linear interpolation
source
VortexLattice.wing_to_surface_panelsFunction
wing_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 section
  • yle: leading edge y-coordinate of each airfoil section
  • zle: leading edge z-coordinate of each airfoil section
  • chord: chord length of each airfoil section
  • theta: twist of each airfoil section
  • phi: dihedral angle of each airfoil section, defined by a right hand rotation about the x-axis
  • ns: number of spanwise panels
  • nc: number of chordwise panels
  • fc: (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 to false
  • 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 to Cosine()
  • spacing_c: chordwise discretization scheme, defaults to Uniform()
  • interp_s: interpolation function between spanwise stations, defaults to linear interpolation
source
VortexLattice.lifting_line_geometryFunction
lifting_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 and ns 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 coordinates
  • c: Vector with length equal to the number of surfaces, with each element being a vector of length ns+1 which contains the chord lengths at each lifting line coordinate.
source
VortexLattice.read_degengeomFunction
`read_degengeom(filename::String)`

Read all geometry components from a DegenGeom file written out by OpenVSP

Arguments

  • filename::String: DegenGeom filename

Returns

  • comp: Vector of vsp.VSPComponent objects
source
VortexLattice.import_vspFunction
`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: Single VSPGeom.VSPComponent object
  • geomType::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 corners
  • surface: Array with dimensions (i, j) containing generated panels
source
VortexLattice.translateFunction
translate(panel::SurfacePanel, r)

Return a copy of panel translated the distance specified by vector r

source
translate(surface, r)

Return a copy of the panels in surface translated the distance specified by vector r

source
translate(surfaces, r)

Return a copy of the surfaces in surfaces translated the distance specified by vector r

source
translate(grid, r)

Return a copy of the grid points in grid translated the distance specified by vector r

source
VortexLattice.translate!Function
translate!(surface, r)

Translate the panels in surface the distance specified by vector r

source
translate!(surfaces, r)

Translate the surfaces in surfaces the distance specified by vector r

source
translate!(grid, r)

Translate the grid points in grid the distance specified by vector r

source
VortexLattice.rotateFunction
rotate(panel::SurfacePanel, R, r = [0,0,0])

Return a copy of panel rotated about point r using the rotation matrix R

source
rotate(surface, R, r = [0,0,0])

Return a copy of the panels in surface rotated about point r using the rotation matrix R

source
rotate(surfaces, R, r = [0,0,0])

Return a copy of the surfaces in surfaces rotated about point r using the rotation matrix R

source
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

source
VortexLattice.rotate!Function
rotate!(surface, R, r = [0,0,0])

Rotate the panels in surface about point r using the rotation matrix R

source
rotate!(surfaces, R, r = [0,0,0])

Rotate the surfaces in surfaces about point r using the rotation matrix R

source
rotate!(grid, R, r = [0,0,0])

Rotate the grid points in grid about point r using the rotation matrix R

source

Reference Parameters and Frames

VortexLattice.ReferenceType
Reference(S, c, b, r)

Reference quantities.

Arguments

  • S: reference area
  • c: reference chord
  • b: reference span
  • r: reference location for all rotations/moments
  • V: reference velocity (magnitude)
source
VortexLattice.StabilityType
Stability <: AbstractFrame

Reference frame rotated from the body frame about the y-axis to be aligned with the freestream alpha.

source
VortexLattice.WindType
Wind <: AbstractFrame

Reference frame rotated to be aligned with the freestream alpha and beta

source

Freestream Parameters

VortexLattice.FreestreamType
Freestream([Vinf,] alpha, beta, Omega)

Defines the freestream and rotational velocity properties.

Arguments

  • Vinf: Freestream velocity
  • alpha: 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)
source
VortexLattice.trajectory_to_freestreamFunction
trajectory_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 step
  • Ydot = zeros(length(dt)): Global frame y-velocity for each time step
  • Zdot = zeros(length(dt)): Global frame z-velocity for each time step
  • p = zeros(length(dt)): Angular velocity about x-axis for each time step
  • q = zeros(length(dt)): Angular velocity about y-axis for each time step
  • r = zeros(length(dt)): Angular velocity about z-axis for each time step
  • phi0 = 0: Roll angle for initial time step
  • theta0 = 0: Pitch angle for initial time step
  • psi0 = 0: Yaw angle for initial time step
source

Performing an Analysis

VortexLattice.SystemType
System{TF}

Contains pre-allocated storage for internal system variables.

Fields:

  • AIC: Aerodynamic influence coefficient matrix from the surface panels
  • w: Normal velocity at the control points from external sources and wakes
  • Γ: Circulation strength of the surface panels
  • V: Velocity at the wake vertices for each surface
  • surfaces: Surfaces, represented by matrices of surface panels
  • properties: Surface panel properties for each surface
  • wakes: Wake panel properties for each surface
  • trefftz: Trefftz panels associated with each surface
  • reference: Pointer to reference parameters associated with the system (see Reference)
  • freestream: Pointer to current freestream parameters associated with the system (see Freestream)
  • symmetric: Flags indicating whether each surface is symmetric across the X-Z plane
  • nwake: Number of chordwise wake panels to use from each wake in wakes,
  • 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 vortices
  • 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.
  • near_field_analysis: Flag indicating whether a near field analysis has been performed for the current system state
  • derivatives: Flag indicating whether the derivatives with respect to the freestream variables have been calculated
  • dw: Derivatives of the R.H.S. with respect to the freestream variables
  • : Derivatives of the circulation strength with respect to the freestream variables
  • dproperties: Derivatives of the panel properties with respect to the freestream variables
  • wake_shedding_locations: Wake shedding locations for each surface
  • Vcp: Velocity due to surface motion at the control points
  • Vh: Velocity due to surface motion at the horizontal bound vortex centers
  • Vv: Velocity due to surface motion at the vertical bound vortex centers
  • Vte: Velocity due to surface motion at the trailing edge vertices
  • dΓdt: Derivative of the circulation strength with respect to non-dimensional time
source
VortexLattice.SystemMethod
System([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 by surface
  • 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 (see SurfacePanel) of shape (nc, ns) where nc is the number of chordwise panels and ns 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.
source
VortexLattice.steady_analysisFunction
steady_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
    or
    • Vector of matrices of shape (nc, ns) containing surface panels (see
    SurfacePanel) where nc is the number of chordwise panels and ns is the number of spanwise panels
  • reference: Reference parameters (see Reference)
  • freestream: Freestream parameters (see Freestream)

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 to false for each surface
  • wakes: Matrix of wake panels (see WakePanel) for each surface. Each matrix has shape (nw, ns) where nw is the number of chordwise wake panels and ns is the number of spanwise panels for each surface, defaults to no wake panels for each surface
  • nwake: Number of chordwise wake panels to use from each wake in wakes, defaults to all wake panels for each surface
  • 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. By default, all surfaces are assigned their own IDs
  • wake_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 to true for each surface.
  • trailing_vortices: Flags to enable/disable trailing vortices for each surface, defaults to true for each surface
  • xhat: 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 to true. 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 to true.
  • derivatives: Flag indicating whether the derivatives with respect to the freestream variables should be calculated. Defaults to true.
source
VortexLattice.unsteady_analysisFunction
unsteady_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, ns) where nc is the number of chordwise panels and ns is the number of spanwise panels. Alternatively, a vector containing surface shapes/positions at each time step (including at t=0) may be provided to model moving/deforming lifting surfaces.
  • reference: Reference parameters (see Reference)
  • freestream: Freestream parameters for each time step (see Freestream)
  • 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 to false for each surface
  • initial_wakes: Vector of initial wakes corresponding to each surface, represented by matrices of wake panels (see WakePanel) of shape (nw, ns) where nw is the number of chordwise wake panels and ns is the number of spanwise panels. Defaults to no wake panels for each surface
  • initial_circulation: Vector containing the initial circulation of all surface panels in the system. Defaults to zeros(N) where N is the total number of surface panels in surfaces.
  • nwake: Maximum number of wake panels in the chordwise direction for each surface. Defaults to length(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 to true for each surface.
  • save: Time indices at which to save the time history, defaults to 1: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 to true. 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 to true.
  • derivatives: Flag indicating whether the derivatives with respect to the freestream variables should be calculated for the final time step. Defaults to true.
source

Near Field Forces and Moments

VortexLattice.PanelPropertiesType
PanelProperties

Panel specific properties calculated during the vortex lattice method analysis.

Fields

  • gamma: Vortex ring circulation strength, normalized by the reference velocity
  • velocity: Local velocity at the panel's bound vortex center, normalized by the reference velocity
  • cfb: Net force on the panel's bound vortex, as calculated using the Kutta-Joukowski theorem, normalized by the reference dynamic pressure and area
  • cfl: 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 area
  • cfr: 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
source
VortexLattice.get_surface_propertiesFunction
get_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

source
VortexLattice.body_forcesMethod
body_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 type System which holds system properties

Keyword Arguments

source
VortexLattice.body_forces_historyFunction
body_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 type System which holds system properties
  • surface_history: Vector of surfaces at each time step, where each surface is represented by a matrix of surface panels (see SurfacePanel) of shape (nc, ns) where nc is the number of chordwise panels and ns is the number of spanwise panels
  • property_history: Vector of surface properties for each surface at each time step, where surface properties are represented by a matrix 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

Keyword Arguments

source
VortexLattice.lifting_line_coefficientsFunction
lifting_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 type System 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 coordinates
  • c: Vector with length equal to the number of surfaces, with each element being a vector of length ns+1 which contains the chord lengths at each lifting line coordinate.

Keyword Arguments

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.
source

Far Field Drag

VortexLattice.far_field_dragFunction
far_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
source
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 (see TrefftzPanel)
  • sending: Vector of sending Trefftz panels (see TrefftzPanel)
  • reference: Reference parameters (see Reference)
  • symmetric: Flag indicating whether a mirror image of the panels in surface, should be used when calculating induced velocities
source

Body and Stability Derivatives

VortexLattice.body_derivativesFunction
body_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 type System which holds system properties
source
VortexLattice.stability_derivativesFunction
stability_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 type System which holds system properties
source

Visualization

VortexLattice.write_vtkFunction
write_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 files
  • surfaces:
    • Vector of grids of shape (3, nc+1, ns+1) which represent lifting surfaces
    or
    • Vector of matrices of shape (nc, ns) containing surface panels (see
    SurfacePanel) where nc is the number of chordwise panels and ns is the number of spanwise panels
  • wakes: (optional) Vector of wakes corresponding to each surface, represented by matrices of wake panels (see WakePanel) of shape (nw, ns) where nw is the number of chordwise wake panels and ns is the number of spanwise panels.
  • surface_properties: (optional) 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

Keyword Arguments:

  • symmetric: (required if surface_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 to true when wake panels are absent, false otherwise
  • xhat: Direction in which trailing vortices extend if used. Defaults to [1, 0, 0].
  • wake_length: Distance to extend trailing vortices. Defaults to 10
  • metadata: Dictionary of metadata to include in generated files
source
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 files
  • surface_history: Vector of surfaces at each time step, where each surface is represented by a matrix of surface panels (see SurfacePanel) of shape (nc, ns) where nc is the number of chordwise panels and ns is the number of spanwise panels
  • property_history: Vector of surface properties for each surface at each time step, where surface properties are represented by a matrix 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
  • wake_history: Vector of wakes corresponding to each surface at each time step, where each wake is represented by a matrix of wake panels (see WakePanel) of shape (nw, ns) where nw is the number of chordwise wake panels and ns is the number of spanwise panels.
  • dt: Time step vector

Keyword Arguments:

  • symmetric: (required if properties 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 10
  • metadata: Dictionary of metadata to include in generated files
source

Private API

Geometry

VortexLattice.linearinterpFunction

linearinterp(eta, rstart, rend)

Linearly interpolate between rstart and rend where eta is the fraction between 0 (rstart) and 1 (rend)

source
VortexLattice.spanwise_spacingFunction
spanwise_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.

source
VortexLattice.chordwise_spacingFunction
chordwise_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.

source
VortexLattice.interpolate_gridFunction
interpolate_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 direction dir (0 <= eta <= 1)
  • interp: Interpolation method of form ypt = f(x,y,xpt)
  • xdir: Independent variable direction, defaults to arc length
  • ydir: Dependent variable direction xyz (i=1, j=2)
source
VortexLattice.repeated_trailing_edge_pointsFunction
repeated_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.

source

Surface Panels

VortexLattice.get_core_sizeFunction
get_core_size(panel::SurfacePanel)

Return the core size (smoothing parameter) corresponding to the vortex ring associated with panel

source
VortexLattice.left_centerFunction
left_center(panel::SurfacePanel)

Return the center of the left bound vortex of the vortex ring associated with panel

source
VortexLattice.right_centerFunction
right_center(panel::SurfacePanel)

Return the center of the right bound vortex of the vortex ring associated with panel

source
VortexLattice.top_vectorFunction
top_vector(panel::SurfacePanel)

Return the path of the top bound vortex of the vortex ring associated with panel

source

Wake Panels

VortexLattice.update_wake_shedding_locations!Function
update_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 (see WakePanel) of shape (nw, ns) where nw is the number of chordwise wake panels and ns 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 (see SurfacePanel of shape (nc, ns) where nc is the number of chordwise panels and ns is the number of spanwise panels
  • reference: Reference parameters (see Reference)
  • freestream: Freestream parameters (see Freestream)
  • dt: Time step (seconds)
  • additional_velocity: Function defining additional velocity field
  • Vte: Velocity experienced at the trailing edge due to surface motion.
  • nwake: Number of chordwise wake panels to use from each wake in wakes
  • eta: Time step fraction used to define separation between trailing edge and wake shedding location. Typical values range from 0.2-0.3.
source
VortexLattice.get_wake_velocities!Function
get_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 in wakes
  • surfaces: Vector of surfaces, represented by matrices of surface panels (see SurfacePanel of shape (nc, ns) where nc is the number of chordwise panels and ns is the number of spanwise panels
  • wakes: Vector of wakes corresponding to each surface, represented by matrices of wake panels (see WakePanel) of shape (nw, ns) where nw is the number of chordwise wake panels and ns is the number of spanwise panels.
  • reference: Reference parameters (see Reference)
  • freestream: Freestream parameters (see Freestream)
  • Γ: Circulation of all surface panels stored in a single vector
  • additional_velocity: Function defining additional velocity field
  • Vte: Velocity at the trailing edge vertices on each surface due to surface motion
  • symmetric: (required) Flag for each surface indicating whether a mirror image across the X-Z plane should be used when calculating induced velocities
  • repeated_points: 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. See repeated_trailing_edge_points
  • nwake: Number of chordwise wake panels to use from each wake in wakes, defaults to all provided wake panels
  • 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 to true 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 to true for each surface
  • xhat: Direction in which to shed trailing vortices, defaults to [1, 0, 0]
source
VortexLattice.translate_wakeFunction
translate_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 type WakePanel)
  • wake_velocities: Matrix containing the velocities at each of the four corners of panel
  • dt: Time step (seconds)
source
VortexLattice.translate_wake!Function
translate_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 (see WakePanel) of shape (nw, ns) where nw is the number of chordwise wake panels and ns is the number of spanwise panels, defaults to no wake panels
  • wake_velocities: Velocities at each of the vertices corresponding to the wake panels in wake
  • dt: Time step

Keyword Arguments

  • nwake: Number of chordwise wake panels to use from wake, defaults to all provided wake panels
source
VortexLattice.shed_wake!Function
shed_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 (see WakePanel) of shape (nw, ns) where nw is the number of chordwise wake panels and ns is the number of spanwise panels
  • wake_shedding_locations: Vector of length ns which stores the coordinates where wake panels are shed from the trailing edge of surface.
  • wake_velocities: Velocities at each of the vertices corresponding to the wake panels in wake
  • dt: Time step (seconds)
  • surface: Matrix of surface panels (see SurfacePanel) of shape (nc, ns) where nc is the number of chordwise panels and ns is the number of spanwise panels
  • Γ: Circulation strength of each surface panel in surface
  • nwake: Number of chordwise wake panels to use from wake, defaults to all provided wake panels
source
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 (see WakePanel) of shape (nw, ns) where nw is the number of chordwise wake panels and ns 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 in wake
  • dt: Time step (seconds)
  • surfaces: Vector of surfaces, represented by matrices of surface panels (see SurfacePanel of shape (nc, ns) where nc is the number of chordwise panels and ns is the number of spanwise panels
  • Γ: Circulation strength of each surface panel in surfaces
  • nwake: Number of chordwise wake panels to use from each wake in wakes, defaults to all provided wake panels
source

Induced Velocity

VortexLattice.bound_induced_velocityFunction
bound_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

source
VortexLattice.trailing_induced_velocityFunction
trailing_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.

source
VortexLattice.ring_induced_velocityFunction
ring_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.

source
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.

source
VortexLattice.influence_coefficients!Function
influence_coefficients!(AIC, surface; kwargs...)

Construct the aerodynamic influence coefficient matrix for a single surface.

Arguments:

  • surface: Matrix of surface panels (see SurfacePanel) of shape (nc, ns) where nc is the number of chordwise panels and ns is the number of spanwise panels

Keyword Arguments

  • symmetric: Flag indicating whether a mirror image of the panels in surface should be used when calculating induced velocities.
  • wake_shedding_locations: Wake shedding locations for the trailing edge panels in surface
  • trailing_vortices: Flag to enable/disable trailing vortices. Defaults to true.
  • xhat: Direction in which trailing vortices are shed if trailing_vortices = true. Defaults to [1, 0, 0]
source
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 (see SurfacePanel of shape (nc, ns) where nc is the number of chordwise panels and ns 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 to false 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 surface
  • wake_shedding_locations: Wake shedding locations for the trailing edge panels of each surface in surfaces
  • trailing_vortices: Flags to indicate whether trailing vortices are used for each surface. Defaults to true for each surface.
  • xhat: Direction in which trailing vortices are shed if trailing_vortices = true. Defaults to [1, 0, 0]
source
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 to true
  • symmetric: Flag indicating whether sending panels should be mirrored across the X-Z plane
  • wake_shedding_locations: Wake shedding locations for the trailing edge panels in sending
  • trailing_vortices: Indicates whether trailing vortices are used. Defaults to true.
  • xhat: Direction in which trailing vortices are shed if trailing_vortices = true. Defaults to [1, 0, 0]
source
VortexLattice.update_trailing_edge_coefficients!Function
update_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) where nc is the number of chordwise panels and ns is the number of spanwise panels

Keyword Arguments

  • symmetric: Flag indicating whether a mirror image of the panels in surface should be used when calculating induced velocities.
  • trailing_vortices: Flag to enable/disable trailing vortices.
  • xhat: Direction in which trailing vortices are shed if trailing_vortices = true. Defaults to [1, 0, 0]
source
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) where nc is the number of chordwise panels and ns 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 to false 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 surface
  • trailing_vortices: Flags to indicate whether trailing vortices are used for each surface. Defaults to true for each surface.
  • xhat: Direction in which trailing vortices are shed if trailing_vortices = true. Defaults to [1, 0, 0]
source
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 to true
  • symmetric: Flag indicating whether sending panels should be mirrored across the X-Z plane
  • wake_shedding_locations: Wake shedding locations for the trailing edge panels in sending
  • trailing_vortices: Indicates whether trailing vortices are used. Defaults to true.
  • xhat: Direction in which trailing vortices are shed if trailing_vortices = true. Defaults to [1, 0, 0]
source
VortexLattice.induced_velocityFunction
induced_velocity(rcp, surface, Γ; kwargs...)

Compute the velocity induced by the grid of panels in surface at control point rcp

source
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.

source
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

source
induced_velocity(rcp, wake::AbstractMatrix{<:WakePanel}; kwargs...)

Compute the velocity induced by the grid of wake panels in wake at control point rcp

source
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

source
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 to size(surface, 1)
  • ns: Number of panels in the spanwise direction. Defaults to size(surface, 2)
  • finite_core: Flag indicating whether the finite core model should be used
  • symmetric: Flag indicating whether sending panels should be mirrored across the X-Z plane
  • wake_shedding_locations: Wake shedding locations for the trailing edge panels in surface
  • trailing_vortices: Indicates whether trailing vortices are used. Defaults to true.
  • xhat: Direction in which trailing vortices are shed if trailing_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 if trailing_vortices = true
  • skip_top: Tuple containing panel indices whose top bound vortex is coincident with rcp and should therefore be skipped.
  • skip_bottom: Tuple containing panel indices whose bottom bound vortex is coincident with rcp and should therefore be skipped.
  • skip_left: Tuple containing panel indices whose left bound vortex is coincident with rcp and should therefore be skipped.
  • skip_right: Tuple containing panel indices whose right bound vortex is coincident with rcp and should therefore be skipped.
  • skip_left_trailing: Tuple containing panel indices whose left trailing vortex is coincident with rcp and should therefore be skipped.
  • skip_right_trailing: Tuple containing panel indices whose right trailing vortex is coincident with rcp and should therefore be skipped.
source
VortexLattice.induced_velocity_derivativesFunction
induced_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

source
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

source

Freestream

VortexLattice.get_surface_velocities!Function
get_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.

source

Circulation

VortexLattice.normal_velocity!Function
normal_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.

source
VortexLattice.normal_velocity_derivatives!Function
normal_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).

source

Time-Domain Analysis

VortexLattice.propagate_system!Function
propagate_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 type system which contains the current system state
  • surfaces: 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 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. See repeated_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 argument surfaces 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.
source

Near-Field Analysis

VortexLattice.near_field_forces!Function
near_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.

source
VortexLattice.near_field_forces_derivatives!Function
near_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.

source
VortexLattice.body_forcesMethod
body_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 (see SurfacePanel) of shape (nc, ns) where nc is the number of chordwise panels and ns is the number of spanwise panels
  • properties: Surface properties for each surface, where surface properties for each surface are represented by a matrix 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
  • reference: Reference parameters (see Reference)
  • 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 velocities
  • frame: frame in which to return CF and CM, options are Body() (default), Stability(), and Wind()
source
VortexLattice.body_forces_derivativesFunction
body_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 type System which holds system properties
source
VortexLattice.body_to_frameFunction
body_to_frame(CF, CM, reference, freestream, frame)

Transform the coefficients CF and CM from the body frame to the frame specified in frame

source

Far-Field

VortexLattice.trefftz_panel_induced_dragFunction
trefftz_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 of sending should be used when calculating the induced drag
source

Visualization

VortexLattice.write_vtk!Function
write_vtk!(vtmfile, surface, [surface_properties]; kwargs...)

Writes geometry to Paraview files for visualization.

Arguments

  • vtmfile: Multiblock file handle
  • surface: Matrix of surface panels (see SurfacePanel) of shape (nc, ns) where nc is the number of chordwise panels and ns is the number of spanwise panels
  • surface_properties: (optional) Matrix of panel properties for each non-wake panel where each element of the matrix is of type PanelProperties.

Keyword Arguments:

  • symmetric: (required if properties 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 vortices
  • xhat = [1, 0, 0]: Direction in which trailing vortices extend if used
  • wake_length = 10: Distance to extend trailing vortices
  • wake_circulation = zeros(size(surfaces, 2)): Contribution to the trailing edge circulation from the wake attached to this surface
  • metadata = Dict(): Dictionary of metadata to include in generated files
source
write_vtk!(vtmfile, wake; kwargs...)

Writes geometry to Paraview files for visualization.

Arguments

  • vtmfile: Paraview file handle
  • wake: Matrix of wake panels (see WakePanel) of shape (nw, ns) where nw is the number of chordwise wake panels and ns 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 vortices
  • xhat = [1, 0, 0]: Direction in which trailing vortices extend if used
  • wake_length = 10: Distance to extend trailing vortices
  • surface_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
source

Index