Introduction
The SystemStructure
provides a flexible framework for defining the physical structure of airborne wind energy (AWE) systems using discrete mass-spring-damper models. This structure can represent many different AWE system configurations, from simple single-line kites to complex multi-wing systems with intricate bridle networks.
The SystemStructure
serves as input to the SymbolicAWEModel
, which is based on ModelingToolkit and automatically generates symbolic differential algebraic equations from the structural definition.
Workflow
- Define system components (
Point
,Segment
,Group
, etc.) - Assemble into a
SystemStructure
- Pass to
SymbolicAWEModel
for automatic MTK model generation - Simulate the resulting symbolic model
Public enumerations
KiteModels.SegmentType
— TypeSegmentType `POWER_LINE` `STEERING_LINE` `BRIDLE`
Type of segment.
Elements
- POWER_LINE: Belongs to a power line
- STEERING_LINE: Belongs to a steering line
- BRIDLE: Belongs to the bridle
KiteModels.DynamicsType
— TypeDynamicsType `DYNAMIC` `QUASI_STATIC` `WING` `STATIC`
Enumeration of the models that are attached to a point.
Elements
- DYNAMIC: Belongs to a dynamic tether model
- QUASI_STATIC: Belongs to a quasi static tether model
- WING: Connected to the rigid wing body
- STATIC: Does not change position
Public constructors
KiteModels.SystemStructure
— MethodSystemStructure(name, set; points=Point[], groups=Group[], segments=Segment[],
pulleys=Pulley[], tethers=Tether[], winches=Winch[],
wings=Wing[], transforms=Transform[])
Constructs a SystemStructure object representing a complete kite system using a discrete mass-spring-damper model.
Components
- Points - See
Point
for discrete mass dynamics - Segments - See
Segment
for elastic spring-damper connections - Groups - See
Group
for wing twist deformation modeling - Wings - See
Wing
for rigid body dynamics - Pulleys - See
Pulley
for length redistribution between segments - Tethers - See
Tether
for segment groups with shared unstretched length - Winches - See
Winch
for ground-based tether length control - Transforms - See
Transform
for initial positioning and orientation
Physical Models
- "ram": 4 deformable wing groups, complex pulley bridle system
- "simple_ram": 4 deformable wing groups, direct bridle connections
Arguments
name::String
: Model identifier. "ram" and "simple_ram" are defined inside KiteModels.jl, provide a different name for a custom model.set::Settings
: Configuration parameters (see KiteUtils.Settings)
Returns
SystemStructure
: Complete system ready for building aSymbolicAWEModel
Examples
# Auto-generate from wing geometry
wing = RamAirWing(set)
sys_struct = SystemStructure(set, wing)
# Manual construction
points = [Point(1, [0,0,0], STATIC), Point(2, [0,0,10], DYNAMIC)]
segments = [Segment(1, (1,2), BRIDLE)]
sys_struct = SystemStructure("custom", set; points, segments)
KiteModels.SystemStructure
— Typestruct SystemStructure
A discrete mass-spring-damper representation of a kite system, where point masses connected by elastic segments model the kite and tether dynamics.
Components
Point
: Point masses representing wing attachment points, dynamic bridle/tether points, and fixed ground anchorsGroup
: Collections of points that move together according to wing deformation (twist and trailing edge deflection)Segment
: Spring-damper elements connecting pointsPulley
: Elements that redistribute line lengths between segmentsTether
: Groups of segments with a common unstretched lengthWinch
: Ground-based winches that control tether lengthsWing
: Rigid wing bodies that serve as reference framesTransform
: Spatial transformations for initial positioning and orientation
See the individual component documentation for detailed mathematical models and governing equations.
KiteModels.SymbolicAWEModel
— MethodSymbolicAWEModel(set::Settings, sys_struct::SystemStructure,
vsm_aeros::Vector{<:BodyAerodynamics}=BodyAerodynamics[],
vsm_solvers::Vector{<:VortexStepMethod.Solver}=VortexStepMethod.Solver[])
Constructs a SymbolicAWEModel that can generate ModelingToolkit equations from the discrete mass-spring-damper representation defined in the SystemStructure
. The aerodynamic models provide forces and moments acting on wing components.
Arguments
set::Settings
: Configuration parameters (see KiteUtils.Settings)sys_struct::SystemStructure
: Physical system definition with points, segments, groups, etc.vsm_aeros::Vector{<:BodyAerodynamics}=BodyAerodynamics[]
: Aerodynamic models for each wingvsm_solvers::Vector{<:VortexStepMethod.Solver}=VortexStepMethod.Solver[]
: VSM solvers for aerodynamic calculations
Returns
SymbolicAWEModel
: Model ready for symbolic equation generation viainit_sim!
Example
# Create wing geometry and aerodynamics
set = se()
wing = RamAirWing(set)
aero = BodyAerodynamics([wing])
solver = Solver(aero; solver_type=NONLIN)
# Create system structure
sys_struct = SystemStructure(set, wing)
# Create symbolic model
model = SymbolicAWEModel(set, sys_struct, [aero], [solver])
KiteModels.SymbolicAWEModel
— MethodSymbolicAWEModel(set::Settings)
Constructs a default SymbolicAWEModel with automatically generated components.
This convenience constructor creates a complete AWE model using default configurations:
- Generates a ram-air wing from settings
- Creates aerodynamic model and VSM solver
- Builds system structure based on the wing geometry
- Assembles everything into a ready-to-use symbolic model
Arguments
set::Settings
: Configuration parameters (see KiteUtils.Settings)
Returns
SymbolicAWEModel
: Model ready for symbolic equation generation viainit_sim!
Example
set = se() # Load default settings
model = SymbolicAWEModel(set)
# Initialize and run simulation
init_sim!(model)
for i in 1:1000
next_step!(model)
end
KiteModels.Point
— MethodPoint(idx, pos_cad, type; wing_idx=1, vel_w=zeros(KVec3), transform_idx=1)
Constructs a Point object. A point can be of four different DynamicsType
s:
STATIC
: the point doesn't move. $\ddot{\mathbf{r}} = \mathbf{0}$DYNAMIC
: the point moves according to Newton's second law. $\ddot{\mathbf{r}} = \mathbf{F}/m$QUASI_STATIC
: the acceleration is constrained to be zero, by solving a nonlinear problem. $\mathbf{F}/m = \mathbf{0}$WING
: the point has a static position in the rigid body wing frame. $\mathbf{r}_w = \mathbf{r}_{wing} + \mathbf{R}_{b\rightarrow w} \mathbf{r}_b$
where:
- $\mathbf{r}$ is the point position vector
- $\mathbf{F}$ is the net force acting on the point
- $m$ is the point mass
- $\mathbf{r}_w$ is the position in world frame
- $\mathbf{r}_{wing}$ is the wing center position
- $\mathbf{R}_{b\rightarrow w}$ is the rotation matrix from body to world frame
- $\mathbf{r}_b$ is the position in body frame
Arguments
idx::Int16
: Unique identifier for the point.pos_cad::KVec3
: Position of the point in the CAD frame.type::DynamicsType
: Dynamics type of the point (STATIC, DYNAMIC, etc.).
Keyword Arguments
wing_idx::Int16=1
: Index of the wing this point is attached to.vel_w::KVec3=zeros(KVec3)
: Initial velocity of the point in world frame.transform_idx::Int16=1
: Index of the transform used for initial positioning.
Returns
Point
: A new Point object.
Example
To create a Point:
point = Point(1, [1.0, 2.0, 3.0], DYNAMIC; wing_idx=1)
KiteModels.Point
— Typemutable struct Point
A point mass.
idx::Int16
transform_idx::Int16
wing_idx::Int16
pos_cad::StaticArraysCore.MVector{3, Float64}
pos_b::StaticArraysCore.MVector{3, Float64}
pos_w::StaticArraysCore.MVector{3, Float64}
vel_w::StaticArraysCore.MVector{3, Float64}
type::DynamicsType
KiteModels.Group
— MethodGroup(idx, point_idxs, vsm_wing::RamAirWing, gamma, type, moment_frac)
Constructs a Group object representing a collection of points on a kite body that share a common twist deformation.
A Group models the local deformation of a kite wing section through twist dynamics. All points within a group undergo the same twist rotation about the chord vector.
The governing equation is:
\[\begin{aligned} \tau = \underbrace{\sum_{i=1}^{4} r_{b,i} \times (\mathbf{F}_{b,i} \cdot \hat{\mathbf{z}})}_{\text{bridles}} + \underbrace{r_a \times (\mathbf{F}_a \cdot \hat{\mathbf{z}})}_{\text{aero}} \end{aligned}\]
where:
- $\tau$ is the total torque about the twist axis
- $r_{b,i}$ is the position vector of bridle point $i$ relative to the twist center
- $\mathbf{F}_{b,i}$ is the force at bridle point $i$
- $\hat{\mathbf{z}}$ is the unit vector along the twist axis (chord direction)
- $r_a$ is the position vector of the aerodynamic center relative to the twist center
- $\mathbf{F}_a$ is the aerodynamic force at the group's aerodynamic center
The group can have two DynamicsType
s:
DYNAMIC
: the group rotates according to Newton's second law: $I\ddot{\theta} = \tau$QUASI_STATIC
: the rotational acceleration is zero: $\tau = 0$
Arguments
idx::Int16
: Unique identifier for the grouppoint_idxs::Vector{Int16}
: Indices of points that move together with this group's twistvsm_wing::RamAirWing
: Wing geometry object used to extract local chord and spanwise vectorsgamma
: Spanwise parameter (typically -1 to 1) defining the group's location along the wingtype::DynamicsType
: Dynamics type (DYNAMIC for time-varying twist, QUASI_STATIC for equilibrium)moment_frac::SimFloat
: Chordwise position (0=leading edge, 1=trailing edge) about which the group rotates
Returns
Group
: A new Group object with twist dynamics capability
Example
Create a group at mid-span with quarter of the wing moment:
group = Group(1, [1, 2, 3], vsm_wing, 0.0, DYNAMIC, 0.25)
KiteModels.Group
— Typestruct Group
Set of bridle lines that share the same twist angle and trailing edge angle.
idx::Int16
point_idxs::Vector{Int16}
le_pos::StaticArraysCore.MVector{3, Float64}
chord::StaticArraysCore.MVector{3, Float64}
y_airf::StaticArraysCore.MVector{3, Float64}
type::DynamicsType
moment_frac::Float64
twist::Float64
twist_vel::Float64
KiteModels.Segment
— MethodSegment(idx, point_idxs, type; l0=zero(SimFloat), compression_frac=0.1)
Constructs a Segment object representing an elastic spring-damper connection between two points.
The segment follows Hooke's law with damping and aerodynamic drag:
Spring-Damper Force:
\[\mathbf{F}_{spring} = \left[k(l - l_0) - c\dot{l}\right]\hat{\mathbf{u}}\]
Aerodynamic Drag:
\[\mathbf{F}_{drag} = \frac{1}{2}\rho C_d A |\mathbf{v}_a| \mathbf{v}_{a,\perp}\]
Total Force:
\[\mathbf{F}_{total} = \mathbf{F}_{spring} + \mathbf{F}_{drag}\]
where:
- $k = \frac{E \pi d^2/4}{l}$ is the axial stiffness
- $l$ is current length, $l_0$ is unstretched length
- $c = \frac{\xi}{c_{spring}} k$ is damping coefficient
- $\hat{\mathbf{u}} = \frac{\mathbf{r}_2 - \mathbf{r}_1}{l}$ is unit vector along segment
- $\dot{l} = (\mathbf{v}_1 - \mathbf{v}_2) \cdot \hat{\mathbf{u}}$ is extension rate
- $\mathbf{v}_{a,\perp}$ is apparent wind velocity perpendicular to segment
Arguments
idx::Int16
: Unique identifier for the segment.point_idxs::Tuple{Int16, Int16}
: Tuple containing the indices of the two points connected by this segment.type::SegmentType
: Type of the segment (POWERLINE, STEERINGLINE, BRIDLE).
Keyword Arguments
l0::SimFloat=zero(SimFloat)
: Unstretched length of the segment. Calculated from point positions if zero.compression_frac::SimFloat=0.1
: Compression fraction of stiffness for compression behavior.
Returns
Segment
: A new Segment object.
Example
To create a Segment:
segment = Segment(1, (1, 2), BRIDLE; l0=10.0)
KiteModels.Segment
— Typemutable struct Segment
A segment from one point index to another point index.
idx::Int16
point_idxs::Tuple{Int16, Int16}
type::SegmentType
l0::Float64
compression_frac::Float64
diameter::Float64
KiteModels.Pulley
— MethodPulley(idx, segment_idxs, type)
Constructs a Pulley object that enforces length redistribution between two segments.
The pulley constraint maintains constant total length while allowing force transmission:
Constraint Equations:
\[l_1 + l_2 = l_{total} = \text{constant}\]
Force Balance:
\[F_{pulley} = F_1 - F_2\]
Dynamics:
\[m\ddot{l}_1 = F_{pulley} = F_1 - F_2\]
where:
- $l_1, l_2$ are the lengths of connected segments
- $F_1, F_2$ are the spring forces in the segments
- $m = \rho_{tether} \pi (d/2)^2 l_{total}$ is the total mass of both segments
- $\dot{l}_1 + \dot{l}_2 = 0$ (velocity constraint)
The pulley can have two DynamicsType
s:
DYNAMIC
: the length redistribution follows Newton's second law: $m\ddot{l}_1 = F_1 - F_2$QUASI_STATIC
: the forces are balanced instantaneously: $F_1 = F_2$
Arguments
idx::Int16
: Unique identifier for the pulley.segment_idxs::Tuple{Int16, Int16}
: Tuple containing the indices of the two segments connected by this pulley.type::DynamicsType
: Dynamics type of the pulley (DYNAMIC or QUASI_STATIC).
Returns
Pulley
: A new Pulley object.
Example
To create a Pulley:
pulley = Pulley(1, (1, 2), DYNAMIC)
KiteModels.Pulley
— Typemutable struct Pulley
A pulley described by two segments with the common point of the segments being the pulley.
idx::Int16
segment_idxs::Tuple{Int16, Int16}
type::DynamicsType
sum_length::Float64
length::Float64
vel::Float64
KiteModels.Tether
— TypeTether(idx, segment_idxs)
Constructs a Tether object representing a flexible line composed of multiple segments.
A tether enforces a shared unstretched length constraint across all its constituent segments:
Length Constraint:
\[\sum_{i \in \text{segments}} l_{0,i} = L\]
Winch Control: The unstretched tether length is controlled by winch acceleration:
\[\ddot L = \alpha(v, F, u)\]
where:
- $L$ is the tether length
- $l_{0,i}$ is the segment unstretched length
- $\alpha(v, F, u)$ is the winch acceleration function depending on model type
Arguments
idx::Int16
: Unique identifier for the tethersegment_idxs::Vector{Int16}
: Indices of segments that form this tether
Returns
Tether
: A new Tether object
Example
Create a tether from segments 1, 2, and 3:
tether = Tether(1, [1, 2, 3])
KiteModels.Winch
— MethodWinch(idx, model, tether_idxs, tether_length; tether_vel=0.0)
Constructs a Winch object that controls tether length through torque or speed regulation.
Tether Length Control:
\[\ddot{L} = \alpha(v, F, u)\]
where:
- $L$ is the tether length
- $v$ is the reel out velocity (tether extension rate)
- $F$ is the tether force
- $u$ is the applied torque or speed setpoint
- $\alpha(v, F, u)$ is the winch acceleration function depending on model type
where the winch acceleration function f_winch
depends on the winch model type:
- Torque-controlled: Direct torque input with motor dynamics
- Speed-controlled: Velocity regulation with internal control loops
For detailed mathematical models of winch dynamics, motor characteristics, and control algorithms, see the WinchModels.jl documentation.
Arguments
idx::Int16
: Unique identifier for the winch.model::AbstractWinchModel
: The winch model (TorqueControlledMachine, AsyncMachine, etc.).tether_idxs::Vector{Int16}
: Vector containing the indices of the tethers connected to this winch.tether_length::SimFloat
: Initial tether length.
Keyword Arguments
tether_vel::SimFloat=0.0
: Initial tether velocity (reel-out rate).
Returns
Winch
: A new Winch object.
Example
To create a Winch:
winch = Winch(1, TorqueControlledMachine(set), [1, 2], 100.0)
KiteModels.Winch
— Typemutable struct Winch
A set of tethers or just one tether connected to a winch.
idx::Int16
model::AbstractWinchModel
tether_idxs::Vector{Int16}
tether_length::Float64
tether_vel::Float64
KiteModels.Wing
— MethodWing(idx, group_idxs, R_b_c, pos_cad; transform_idx=1, angular_vel=zeros(KVec3),
pos_w=zeros(KVec3), vel_w=zeros(KVec3))
Constructs a Wing object representing a rigid body that serves as a reference frame for attached points and groups.
A Wing provides a rigid body coordinate system for kite components. Points with type == WING
move rigidly with the wing body according to the wing's orientation matrix and position. Groups attached to the wing undergo local deformation (twist) relative to the rigid wing body frame.
Rigid Body Dynamics: The wing follows standard rigid body equations of motion:
\[\begin{aligned} \frac{\delta \mathbf{q}_b^w}{\delta t} &= \frac{1}{2} \Omega(\boldsymbol{\omega}_b) \mathbf{q}_b^w \\ \boldsymbol{\tau}_w &= \mathbf{I} \frac{\delta \boldsymbol{\omega}}{\delta t} + \boldsymbol{\omega}_b \times (\mathbf{I}\boldsymbol{\omega}_b) \end{aligned}\]
where:
- $\mathbf{q}_b^w$ is the quaternion from world to body frame
- $\boldsymbol{\omega}_b$ is the angular velocity in body frame
- $\Omega(\boldsymbol{\omega}_b)$ is the quaternion multiplication matrix
- $\mathbf{I}$ is the inertia tensor in body frame
- $\boldsymbol{\tau}_w$ is the total applied torque to the rigid wing body (aerodynamic + tether forces)
Coordinate Transformations: Points attached to the wing transform as:
\[\mathbf{r}_w = \mathbf{r}_{w} + \mathbf{R}_{b \rightarrow w} \mathbf{r}_b\]
where:
- $\mathbf{r}_w$ is the position in world frame
- $\mathbf{r}_{w}$ is the wing position in world frame
- $\mathbf{R}_{b \rightarrow w}$ is the rotation from body to world frame
- $\mathbf{r}_b$ is the point position in body frame
Arguments
idx::Int16
: Unique identifier for the winggroup_idxs::Vector{Int16}
: Indices of groups attached to this wing that can deform relative to the bodyR_b_c::Matrix{SimFloat}
: Rotation matrix from body frame to CAD frame (3×3 orthogonal matrix)pos_cad::KVec3
: Position of wing center of mass in CAD frame
Keyword Arguments
transform_idx::Int16=1
: Transform used for initial positioning and orientationangular_vel::KVec3=zeros(KVec3)
: Initial angular velocity of the wing in world framepos_w::KVec3=zeros(KVec3)
: Initial position of wing center of mass in world framevel_w::KVec3=zeros(KVec3)
: Initial velocity of wing center of mass in world frame
Special Properties
The wing orientation can be accessed as a quaternion:
quat = wing.orient # Returns quaternion representation of R_b_c
Returns
Wing
: A new Wing object providing a rigid body reference frame
Example
Create a wing with identity orientation and two attached groups:
R_b_c = I(3) # identity matrix
pos_cad = [0.0, 0.0, 0.0]
wing = Wing(1, [1, 2], R_b_c, pos_cad)
KiteModels.Wing
— Typestruct Wing
A rigid wing body that can have multiple groups of points attached to it.
Fields
idx::Int16
: Unique identifier for the winggroup_idxs::Vector{Int16}
: Indices of groups attached to this wingtransform_idx::Int16
: Transform used for initial positioning and orientationR_b_c::Matrix{SimFloat}
: Rotation matrix from body frame to CAD frameangular_vel::KVec3
: Angular velocity of the wing in world framepos_w::KVec3
: Position of wing center of mass in world framepos_cad::KVec3
: Position of wing center of mass in CAD framevel_w::KVec3
: Velocity of wing center of mass in world frame
The wing provides a rigid body reference frame for attached points and groups. Points with type == WING
move rigidly with the wing body according to the wing's orientation matrix R_b_c
and position pos_w
.
Extended help
The wing's orientation can be accessed as a quaternion through the orient
property:
wing = Wing(1, [1,2], I(3), zeros(3))
quat = wing.orient # Returns quaternion representation of R_b_c
KiteModels.Transform
— MethodTransform(idx, elevation, azimuth, heading;
base_point_idx=nothing, base_pos=nothing, base_transform_idx=nothing,
wing_idx=nothing, rot_point_idx=nothing)
Constructs a Transform object that orients system components using spherical coordinates.
All points and wings with matching transform_idx
are transformed together as a rigid body:
- Translation: Translate such that base is at specified base pos
- Rotation 1: Rotate so target is at (elevation, azimuth) relative to base
- Rotation 2: Rotate all components by
heading
around the base-target vector
\[\mathbf{r}_{transformed} = \mathbf{r}_{base} + \mathbf{R}_{heading} \circ \mathbf{R}_{elevation,azimuth}(\mathbf{r} - \mathbf{r}_{base})\]
Arguments
idx::Int16
: Unique identifier for the transformelevation::SimFloat
: Target elevation angle from base (radians)azimuth::SimFloat
: Target azimuth angle from base (radians)heading::SimFloat
: Rotation around base-target vector (radians)
Keyword Arguments
Base Reference (choose one):
base_pos + base_point_idx
: Fixed position and reference pointbase_transform_idx
: Chain to another transform's position
Target Object (choose one):
wing_idx
: Wing to position at (elevation, azimuth)rot_point_idx
: Point to position at (elevation, azimuth)
Returns
Transform
: Transform affecting all components with matchingtransform_idx
Examples
# Position wing and all associated points at 45° elevation, 30° azimuth
transform = Transform(1, deg2rad(45), deg2rad(30), 0.0;
base_pos=[0,0,0], base_point_idx=1, wing_idx=1)
# Chain transforms for multi-kite systems
transform2 = Transform(2, deg2rad(30), deg2rad(45), deg2rad(10);
base_transform_idx=1, wing_idx=2)
KiteModels.Transform
— Typemutable struct Transform
Describes the spatial transformation (position and orientation) of system components relative to a base reference point.
idx::Int16
wing_idx::Union{Nothing, Int16}
rot_point_idx::Union{Nothing, Int16}
base_point_idx::Union{Nothing, Int16}
base_transform_idx::Union{Nothing, Int16}
elevation::Float64
azimuth::Float64
heading::Float64
base_pos::Union{Nothing, StaticArraysCore.MVector{3, Float64}}
Private functions
KiteModels.wing_eqs!
— Functionwing_eqs!(s, eqs, defaults; kwargs...)
Generate the differential equations for wing dynamics including quaternion kinematics, angular velocities and accelerations, and forces/moments.
Arguments
s::SymbolicAWEModel
: The wing system stateeqs
: Current system equationsdefaults
: Default values for variablestether_wing_force
: Forces from tethers on wingtether_wing_moment
: Moments from tethers on wingaero_force_b
: Aerodynamic forces in body frameaero_moment_b
: Aerodynamic moments in body frameω_b
: Angular velocity in body frameR_b_w
: Body to world rotation matrixwing_pos
: Kite position vectorwing_vel
: Kite velocity vectorwing_acc
: Kite acceleration vectorstabilize
: Whether in stabilize mode
Returns
Tuple of updated equations and defaults
KiteModels.reinit!
— Functionreinit!(s::SymbolicAWEModel, solver; prn=true, precompile=false) -> Nothing
Reinitialize an existing kite power system model with new state values. The new state is coming from the init section of the settings, stored in the struct s.set
.
This function performs the following operations:
- If no integrator exists yet:
- Loads a serialized ODEProblem from disk
- Initializes a new ODE integrator
- Generates getter/setter functions for the system
- Converts initial settings to quaternion orientation
- Initializes the point mass system with new positions
- Sets initial values for all state variables
- Reinitializes the ODE integrator
- Updates the linearized aerodynamic model
This is more efficient than init!
as it reuses the existing model structure and only updates the state variables to match the current initial settings.
Arguments
s::SymbolicAWEModel
: The kite power system state objectsolver
: The solver to be usedprn::Bool=true
: Whether to print progress information
Returns
Nothing
Throws
ArgumentError
: If no serialized problem exists (runinit_sim!
first)
KiteModels.scalar_eqs!
— Functionscalar_eqs!(s, eqs; R_b_w, wind_vec_gnd, va_wing_b, wing_pos, wing_vel, wing_acc, twist_angle, twist_ω)
Generate equations for scalar quantities like elevation, azimuth, heading and course angles.
# Arguments
- `s::SymbolicAWEModel`: The wing system state
- `eqs`: Current system equations
- `R_b_w`: Body to world rotation matrix
- `wind_vec_gnd`: Ground wind vector
- `va_wing_b`: Apparent wind velocity in body frame
- `wing_pos`: Kite position vector
- `wing_vel`: Kite velocity vector
- `wing_acc`: Kite acceleration vector
# Returns
- Updated system equations including:
- Heading angle from x-axis
- Elevation angle
- Azimuth angle
- Course angle
- Angular velocities and accelerations
KiteModels.linear_vsm_eqs!
— Functionlinearvsmeqs!(s, eqs; aeroforceb, aeromomentb, groupaeromoment, initvab, twistangle, vawingb, ωb)
Generate linearized aerodynamic equations using the Vortex Step Method (VSM). Uses linearization around current operating point to approximate aerodynamic forces and moments. The Jacobian is computed using the VSM solver.
Arguments
s::SymbolicAWEModel
: The wing system stateeqs
: Current system equationsaero_force_b
: Aerodynamic forces in body frameaero_moment_b
: Aerodynamic moments in body framegroup_aero_moment
: Aerodynamic moments per groupinit_va_b
: Initial apparent wind velocitytwist_angle
: Twist angles per groupva_wing_b
: Apparent wind velocity in body frameω_b
: Angular velocity in body frame
Returns
- Updated system equations including linearized aerodynamics:
- Force and moment calculations
- Group moment distributions
- Jacobian matrix for state derivatives
KiteModels.force_eqs!
— Functionforce_eqs!(s, system, eqs, defaults, guesses; kwargs...)
Generate the force equations for the wing system including spring forces, drag forces, pulley dynamics and winch forces.
Arguments
s::SymbolicAWEModel
: The wing system statesystem::SystemStructure
: The point mass representationeqs
: Current system equationsdefaults
: Default values for variablesguesses
: Initial guesses for variablesR_b_w
: Body to world rotation matrixwing_pos
: Kite position vectorwing_vel
: Kite velocity vectorwind_vec_gnd
: Ground wind vectorgroup_aero_moment
: Aerodynamic moments per grouptwist_angle
: Twist angles per groupstabilize
: Whether in stabilize mode
Returns
Tuple containing:
- Updated equations
- Updated defaults
- Updated guesses
- Tether forces on wing
- Tether moments on wing