Enumerations
VortexStepMethod.Model
— TypeModel VSM
LLT
Enumeration of the implemented model types.
Elements
- VSM: Vortex Step Method
- LLT: Lifting Line Theory
VortexStepMethod.WingType
— TypeWingType `RECTANGULAR` `CURVED` `ELLIPTICAL`
Enumeration of the implemented wing types.
Elements:
- RECTANGULAR
- CURVED
- ELLIPTICAL
VortexStepMethod.AeroModel
— TypeAeroModel LEI_AIRFOIL_BREUKELS
POLAR_VECTORS
POLAR_MATRICES
INVISCID
Enumeration of the implemented aerodynamic models. See also: AeroData
Elements
LEI_AIRFOIL_BREUKELS
: Polynom approximation for leading edge inflatable kitesPOLAR_VECTORS
: Polar vectors as function of alpha (lookup tables with interpolation)POLAR_MATRICES
: Polar matrices as function of alpha and delta (lookup tables with interpolation)- INVISCID
where alpha
is the angle of attack, delta
is trailing edge angle.
VortexStepMethod.PanelDistribution
— TypePanelDistribution LINEAR
COSINE
COSINE_VAN_GARREL
SPLIT_PROVIDED
UNCHANGED
Enumeration of the implemented panel distributions.
Elements
- LINEAR # Linear distribution
- COSINE # Cosine distribution
COSINE_VAN_GARREL
# van Garrel cosine distributionSPLIT_PROVIDED
# Split provided sections- UNCHANGED # Keep original sections
VortexStepMethod.InitialGammaDistribution
— TypeInitialGammaDistribution ELLIPTIC ZEROS
Enumeration of the implemented initial gamma distributions.
Elements
- ELLIPTIC
- ZEROS
VortexStepMethod.SolverStatus
— TypeSolverStatus FEASIBLE INFEASIBLE FAILURE
Enumeration to report back the validity of the result of the solve! function. Used in the VSMSolution struct.
Elements
- FEASIBLE: The gamma distribution is physically feasible
- INFEASIBLE: The gamma distribution is physically infeasible
- FAILURE: The result did not converge within the maximal number of iterations
Basic Vectors
VortexStepMethod.MVec3
— Typeconst MVec3 = MVector{3, Float64}
Basic 3-dimensional vector, stack allocated, mutable.
VortexStepMethod.PosVector
— Typeconst PosVector=Union{MVec3, Vector}
Position vector, either a MVec3
or a Vector
for use in function signatures.
VortexStepMethod.VelVector
— Typeconst VelVector=Union{MVec3, Vector}
Velocity vector, either a MVec3
or a Vector
for use in function signatures.
Aerodynamic data
VortexStepMethod.AeroData
— TypeAeroData= Union{
Nothing,
NTuple{2, Float64},
Tuple{Vector{Float64}, Vector{Float64}, Vector{Float64}, Vector{Float64}},
Tuple{Vector{Float64}, Vector{Float64}, Matrix{Float64}, Matrix{Float64}, Matrix{Float64}}
}
Union of different definitions of the aerodynamic properties of a wing section. See also: AeroModel
- nothing for INVISCID
- (
tube_diameter
, camber) forLEI_AIRFOIL_BREUKELS
- (
alpha_range
,cl_vector
,cd_vector
,cm_vector
) forPOLAR_VECTORS
- (
alpha_range
,delta_range
,cl_matrix
,cd_matrix
,cm_matrix
) forPOLAR_MATRICES
where alpha
is the angle of attack [rad], delta
is trailing edge angle [rad], cl
the lift coefficient, cd
the drag coefficient and cm
the pitching moment coefficient. The camber of a kite refers to the curvature of its airfoil shape. The camber is typically measured as the maximum distance between the mean camber line (the line equidistant from the upper and lower surfaces) and the chord line of the airfoil.
Wing Geometry, Panel and Aerodynamics
A body is constructed of one or more abstract wings. An abstract wing can be a Wing or a RamAirWing. A Wing/ RamAirWing has one or more sections.
VortexStepMethod.Section
— Type@with_kw mutable struct Section
Represents a wing section with leading edge, trailing edge, and aerodynamic properties.
Fields
VortexStepMethod.Section
— MethodSection(LE_point::PosVector, TE_point::PosVector, aero_model)
Create a new wing section with the specified leading edge point, trailing edge point, and aerodynamic model.
Arguments
LE_point::PosVector
: Leading edge point coordinatesTE_point::PosVector
: Trailing edge point coordinatesaero_model::AeroModel
: Aerodynamic model type (e.g., INVISCID, POLAR_VECTORS)
Returns
Section
: A new section with the specified parameters and no aerodynamic data
VortexStepMethod.Wing
— TypeWing
Represents a wing composed of multiple sections with aerodynamic properties.
Fields
n_panels::Int16
: Number of panels in aerodynamic meshn_groups::Int16
: Number of panel groupsspanwise_distribution
::PanelDistribution: PanelDistributionspanwise_direction::MVec3
: Wing span direction vectorsections::Vector{Section}
: Vector of wing sections, see: Sectionrefined_sections::Vector{Section}
: Vector of refined wing sections, see: Sectionremove_nan::Bool
: Wether to remove the NaNs from interpolations or not
VortexStepMethod.Wing
— MethodWing(n_panels::Int;
n_groups=n_panels,
spanwise_distribution::PanelDistribution=LINEAR,
spanwise_direction::PosVector=MVec3([0.0, 1.0, 0.0]),
remove_nan::Bool=true)
Constructor for a Wing struct with default values that initializes the sections and refined sections as empty arrays.
Parameters
n_panels::Int
: Number of panels in aerodynamic meshn_groups::Int
: Number of panel groups in aerodynamic meshspanwise_distribution
::PanelDistribution = LINEAR: PanelDistributionspanwise_direction::MVec3
= MVec3([0.0, 1.0, 0.0]): Wing span direction vectorremove_nan::Bool
: Wether to remove the NaNs from interpolations or not
VortexStepMethod.RamAirWing
— TypeRamAirWing <: AbstractWing
A ram-air wing model that represents a curved parafoil with deformable aerodynamic surfaces.
Core Features
- Curved wing geometry derived from 3D mesh (.obj file)
- Aerodynamic properties based on 2D airfoil data (.dat file)
- Support for control inputs (twist angles and trailing edge deflections)
- Inertial and geometric properties calculation
Notable Fields
n_panels::Int16
: Number of panels in aerodynamic meshn_groups::Int16
: Number of control groups for distributed deformationmass::Float64
: Total wing mass in kggamma_tip::Float64
: Angular extent from center to wing tipinertia_tensor::Matrix{Float64}
: Full 3x3 inertia tensor in the kite body frameT_cad_body::MVec3
: Translation vector from CAD frame to body frameR_cad_body::MMat3
: Rotation matrix from CAD frame to body frameradius::Float64
: Wing curvature radiustheta_dist::Vector{Float64}
: Panel twist angle distributiondelta_dist::Vector{Float64}
: Trailing edge deflection distribution
See constructor RamAirWing(obj_path, dat_path; kwargs...)
for usage details.
VortexStepMethod.RamAirWing
— MethodRamAirWing(obj_path, dat_path; kwargs...)
Create a ram-air wing model from 3D geometry and airfoil data files.
This constructor builds a complete aerodynamic model by:
- Loading or generating wing geometry from the .obj file
- Creating aerodynamic polars from the airfoil .dat file
- Computing inertial properties and coordinate transformations
- Setting up control surfaces and panel distribution
Arguments
obj_path
: Path to .obj file containing 3D wing geometrydat_path
: Path to .dat file containing 2D airfoil profile
Keyword Arguments
crease_frac=0.9
: Normalized trailing edge hinge location (0-1)wind_vel=10.0
: Reference wind velocity for XFoil analysis (m/s)mass=1.0
: Wing mass (kg)n_panels=56
: Number of aerodynamic panels across wingspann_groups=4
: Number of control groups for deformationn_sections=n_panels+1
: Number of spanwise cross-sectionsalign_to_principal=false
: Align body frame to principal axes of inertiaspanwise_distribution=UNCHANGED
: Panel distribution typeremove_nan=true
: Interpolate NaN values in aerodynamic dataalpha_range=deg2rad.(-5:1:20)
: Angle of attack range for polars (rad)delta_range=deg2rad.(-5:1:20)
: Trailing edge deflection range for polars (rad)- prn=true: if info messages shall be printed
Returns
A fully initialized RamAirWing
instance ready for aerodynamic simulation.
Example
# Create a ram-air wing from geometry files
wing = RamAirWing(
"path/to/wing.obj",
"path/to/airfoil.dat";
mass=1.5,
n_panels=40,
n_groups=4
)
VortexStepMethod.BodyAerodynamics
— Type@with_kw mutable struct BodyAerodynamics{P}
Main structure for calculating aerodynamic properties of bodies. Use the constructor to initialize.
Fields
- panels::Vector{Panel}: Vector of Panel structs
- wings::Union{Vector{Wing}, Vector{RamAirWing}}: A vector of wings; a body can have multiple wings
va::MVec3
= zeros(MVec3): A vector of the apparent wind speed, see: MVec3omega
::MVec3 = zeros(MVec3): A vector of the turn rates around the kite body axesgamma_distribution
=zeros(Float64, P): A vector of the circulation of the velocity field; Length: Number of segments. [m²/s]alpha_uncorrected
=zeros(Float64, P): angles of attack per panelalpha_corrected
=zeros(Float64, P): corrected angles of attack per panelstall_angle_list
=zeros(Float64, P): stall angle per panelalpha_array::MVector{P, Float64}
= zeros(Float64, P)v_a_array::MVector{P, Float64}
= zeros(Float64, P)work_vectors
::NTuple{10, MVec3} = ntuple(_ -> zeros(MVec3), 10)AIC::Array{Float64, 3}
= zeros(3, P, P)projected_area::Float64
= 1.0: The area projected onto the xy-plane of the kite body reference frame [m²]y::MVector{P, Float64}
= zeros(MVector{P, Float64})cache::Vector{PreallocationTools.LazyBufferCache{typeof(identity), typeof(identity)}}
= [LazyBufferCache() for _ in 1:5]
The Solver and its results
VortexStepMethod.Solver
— TypeSolver
Main solver structure for the Vortex Step Method.See also: solve
Attributes
General settings
aerodynamic_model_type
::Model = VSM: The model type, see: Model- density::Float64 = 1.225: Air density [kg/m³]
max_iterations
::Int64 = 1500rtol
::Float64 = 1e-5: relative errortol_reference_error
::Float64 = 0.001relaxation_factor
::Float64 = 0.03: Relaxation factor for convergence
Damping settings
is_with_artificial_damping
::Bool = false: Whether to apply artificial dampingartificial_damping
::NamedTuple{(:k2, :k4), Tuple{Float64, Float64}} = (k2=0.1, k4=0.0): Artificial damping parameters
Additional settings
type_initial_gamma_distribution
::InitialGammaDistribution = ELLIPTIC: see: InitialGammaDistributioncore_radius_fraction
::Float64 = 1e-20:- mu::Float64 = 1.81e-5: Dynamic viscosity [N·s/m²]
is_only_f_and_gamma_output
::Bool = false: Whether to only output f and gamma
Solution
sol::VSMSolution = VSMSolution(): The result of calling solve!
VortexStepMethod.VSMSolution
— TypeVSMSolution
Struct for storing the solution of the solve! function. Must contain all info needed by KiteModels.jl
.
Attributes
panel_width_array
::Vector{Float64}: Width of the panels [m]alpha_array
::Vector{Float64}: Angle of attack of each panel relative to the apparent wind [rad]- cl_array::Vector{Float64}: Lift coefficients of the panels [-]
- cd_array::Vector{Float64}: Drag coefficients of the panels [-]
- cm_array::Vector{Float64}: Pitching moment coefficients of the panels [-]
- panel_lift::Vector{Float64}: Lift force of the panels [N]
- panel_drag::Vector{Float64}: Drag force of the panels [N]
- panel_moment::Vector{Float64}: Pitching moment around the spanwise vector of the panels [Nm]
f_body_3D
::Matrix{Float64}: Matrix of the aerodynamic forces (x, y, z vectors) [N]m_body_3D
::Matrix{Float64}: Matrix of the aerodynamic moments [Nm]gamma_distribution
::Union{Nothing, Vector{Float64}}: Vector containing the panel circulations.- force::MVec3: Aerodynamic force vector in KB reference frame [N]
- moment::MVec3: Aerodynamic moments [Mx, My, Mz] around the reference point [Nm]
- force_coeffs::MVec3: Aerodynamic force coefficients [CFx, CFy, CFz] [-]
moment_coeffs
::MVec3: Aerodynamic moment coefficients [CMx, CMy, CMz] [-]moment_dist
::Vector{Float64}: Pitching moments around the spanwise vector of each panel. [Nm]moment_coeff_dist
::Vector{Float64}: Pitching moment coefficient around the spanwise vector of each panel. [-]solver_status
::SolverStatus: enum, see SolverStatus