@meta
CurrentModule = SymbolicAWEModels
API Reference
This page provides a detailed reference for all public functions exported by the SymbolicAWEModels.jl
package.
High-Level Simulation Functions
These functions provide convenient wrappers for running common simulation scenarios.
SymbolicAWEModels.sim!
— Functionsim!(sam, set_values; dt, total_time, vsm_interval, prn)
Run a generic simulation for a given AWE model and a matrix of control inputs.
This is the core simulation function that steps through time, applies control inputs, and logs the system state.
Arguments
sam::SymbolicAWEModel
: Initialized AWE model.set_values::Matrix{Float64}
: A matrix of external control torques [Nm] to be applied at each time step. The number of rows must equal the number of simulation steps, and the number of columns must equal the number of winches.
Keywords
dt::Float64
: Time step [s]. Defaults to1/sam.set.sample_freq
.total_time::Float64
: Total simulation duration [s]. Defaults to 10.0.vsm_interval::Int
: Interval for the value state machine updates. Defaults to 3.prn::Bool
: If true, prints a performance summary upon completion. Defaults to true.
Details
The set_values
represent external control torques. Internally, the function calculates the total torque for each winch by adding the external torque to the torque generated by the existing tether force on the drum. This total torque is then used as the set-point for the winch controller in the subsequent step.
Returns
SysLog
: A log object containing the time-series data of the simulation.
SymbolicAWEModels.sim_oscillate!
— Functionsim_oscillate!(sam; dt, total_time, steering_freq, steering_magnitude, vsm_interval, bias, prn)
Run a simulation with oscillating steering input on the given AWE model.
This function is a high-level wrapper around sim!
that generates sinusoidal control inputs to simulate figure-eight maneuvers or steering tests.
Arguments
sam::SymbolicAWEModel
: Initialized AWE model.
Keywords
dt::Float64
: Time step [s]. Defaults to1/sam.set.sample_freq
.total_time::Float64
: Simulation duration [s]. Defaults to 10.0.steering_freq::Float64
: Steering oscillation frequency [Hz]. Defaults to 0.5.steering_magnitude::Float64
: Steering torque amplitude [Nm]. This torque is applied positively to one steering winch and negatively to the other. Defaults to 10.0.vsm_interval::Int
: Value state machine interval. Defaults to 3.bias::Float64
: Phase bias [rad] for the steering cosine wave. Defaults to 0.13.prn::Bool
: If true, prints progress and summary information. Defaults to false.
Returns
SysLog
: A log object containing the time-series data of the simulation.
SymbolicAWEModels.sim_turn!
— Functionsim_turn!(sam; dt, total_time, steering_time, steering_magnitude, vsm_interval, prn)
Run a simulation with a constant steering input for a specified duration.
This function is a wrapper around sim!
that generates a step input for steering, holds it for a duration of steering_time
, and then returns to zero. It is useful for analyzing the step response of the kite's turning behavior.
Arguments
sam::SymbolicAWEModel
: Initialized AWE model.
Keywords
dt::Float64
: Time step [s]. Defaults to1/sam.set.sample_freq
.total_time::Float64
: Total simulation duration [s]. Defaults to 10.0.steering_time::Float64
: Duration of the constant steering input [s]. Defaults to 2.0.steering_magnitude::Float64
: Steering torque amplitude [Nm]. Applied positively to one steering winch and negatively to the other. Defaults to 10.0.vsm_interval::Int
: Value state machine interval. Defaults to 3.prn::Bool
: If true, prints progress and summary information. Defaults to false.
Returns
SysLog
: A log object containing the time-series data of the simulation.
Low-Level Simulation and Analysis
These functions provide direct control over the simulation and tools for model analysis.
KiteUtils.init!
— Functioninit!(s::SymbolicAWEModel; solver, adaptive, prn, precompile, remake, reload, lin_outputs) -> ODEIntegrator
Initialize a kite power system model, creating or loading a compiled version.
If a serialized model exists for the current configuration, it will load that model (fast path). Otherwise, it will create a new model from scratch (slow path).
Arguments
s::SymbolicAWEModel
: The kite system state object.
Keyword Arguments
solver
: Solver algorithm to use. Ifnothing
, defaults based ons.set.solver
.adaptive::Bool=true
: Use adaptive time stepping.prn::Bool=false
: Print progress information.precompile::Bool=false
: Build a generic problem for precompilation.remake::Bool=false
: Force the system to be rebuilt, even if a serialized model exists.reload::Bool=false
: Force the system to reload the serialized model from disk.lin_outputs::Vector{Num}=nothing
: List of symbolic variables for linearization.
Returns
OrdinaryDiffEqCore.ODEIntegrator
: The initialized ODE integrator.
KiteUtils.next_step!
— Functionnext_step!(s::SymbolicAWEModel, integrator::ODEIntegrator; set_values, dt, vsm_interval)
Take a simulation step, using the provided integrator.
next_step!(s::SymbolicAWEModel; set_values, dt, vsm_interval)
Take a simulation step forward in time.
Arguments
s::SymbolicAWEModel
: The kite power system state object.
Keyword Arguments
set_values=nothing
: New values for the control inputs. Ifnothing
, the current values are used.dt=1/s.set.sample_freq
: Time step size [s].vsm_interval=1
: Interval (in steps) to re-linearize the VSM model. If 0, it is not re-linearized.
SymbolicAWEModels.find_steady_state!
— Functionfind_steady_state!(s::SymbolicAWEModel, integ=s.integrator; t=1.0, dt=1/s.set.sample_freq)
Run the simulation for a short period to allow the system to settle.
During this period, the winches are braked and the wing's elevation and azimuth angles are fixed, but it is free to move radially (in distance). This allows the dynamic components of the bridle and tethers to settle into a stable, steady-state equilibrium before starting a maneuver or analysis.
Arguments
s::SymbolicAWEModel
: The model to be stabilized.integ
: The integrator to use. Defaults tos.integrator
.
Keywords
t::Float64=1.0
: The duration [s] for which to run the settling simulation.dt::Float64
: The time step [s] for the settling simulation.
SymbolicAWEModels.linearize!
— Functionlinearize!(s::SymbolicAWEModel; set_values=s.get_set_values(s.integrator)) -> LinType
Compute the full state-space linearization of the model around the current operating point.
This function uses the LinearizationProblem
generated by ModelingToolkit.jl
to calculate the A, B, C, and D matrices for the complete, high-order system.
Arguments
s::SymbolicAWEModel
: The model to linearize.
Keywords
set_values
: The control input vectoru
around which to linearize.
Returns
LinType
: A NamedTuple(A, B, C, D)
containing the state-space matrices.
SymbolicAWEModels.copy_to_simple!
— Functioncopy_to_simple!(sam, tether_sam, simple_sam; prn=true)
Simplify a detailed AWE model into a 1-segment tether model.
This high-level function orchestrates the simplification process:
- Calculates the equivalent axial stiffness and damping of the detailed model (
tether_sam
) by analyzing its step response. - Assigns these calculated properties to the single-segment tethers of the simple model (
simple_sam
). - Copies the dynamic state (wing position, orientation, tether attachment points, etc.) from the detailed model (
sam
) to the simple model (simple_sam
). - Reinitializes the simple model's integrator to apply the new state.
Arguments
sam::SymbolicAWEModel
: The detailed source model, used as a reference for state.tether_sam::SymbolicAWEModel
: A copy of the detailed model, used to perform the step response test.simple_sam::SymbolicAWEModel
: The destination simple model to be updated.
Keywords
prn::Bool=true
: If true, enables printing during the process.
Returns
SymbolicAWEModel
: The updated simple modelsimple_sam
.
copy_to_simple!(sys::SystemStructure, ssys::SystemStructure)
Copy the dynamic state from a detailed SystemStructure
to a simplified one.
This is a low-level utility that maps the state of a complex model (e.g., "ram" with 4 groups and bridle pulleys) to a simpler model (e.g., "simple_ram" with 2 groups and direct connections). It ensures the simplified model matches the key dynamic properties of the detailed one.
Arguments
sys::SystemStructure
: The sourceram
model structure.ssys::SystemStructure
: The destinationsimple_ram
model structure.
SymbolicAWEModels.simple_linearize!
— Functionsimple_linearize!(s::SymbolicAWEModel; tstab=10.0) -> LinType
Compute a simplified, low-order state-space model by numerically linearizing the full simulation.
This function performs system identification by perturbing the states and inputs of the full nonlinear model and observing the effect on the state derivatives and outputs. It runs the simulation for a short duration (tstab
) after each perturbation to find the steady-state response.
Arguments
s::SymbolicAWEModel
: The model to linearize.
Keywords
tstab::Float64=10.0
: The simulation time [s] to run after each perturbation to reach a steady state.
Returns
LinType
: A NamedTuple(A, B, C, D)
containing the identified low-order state-space matrices.
System Structure Constructors
These functions are used to procedurally generate predefined SystemStructure
topologies.
SymbolicAWEModels.create_ram_sys_struct
— Functioncreate_ram_sys_struct(set::Settings)
Create a SystemStructure
for the primary "ram" model with a stability-enhancing bridle.
This is the main, detailed model configuration. It differs from the 4_attach
version by having each of its four Group
sections defined by three deforming bridle points and one statically attached (non-deforming) point. This design improves stability by ensuring the kite's z-axis remains aligned with the bridle system.
The model features:
- A flexible wing with 4 deformable groups (3 deforming points + 1 static point each).
- A complex bridle system with pulleys.
- Four main tethers and three winches.
Arguments
set::Settings
: Configuration parameters defining the kite's geometry and properties.
Returns
SystemStructure
: A newSystemStructure
object representing the "ram" model.
SymbolicAWEModels.create_tether_sys_struct
— Functioncreate_tether_sys_struct(set::Settings; axial_stiffness, axial_damping)
Create a simplified SystemStructure
for testing tether dynamics.
This model consists of only four independent tethers, each represented by a dynamic point mass connected to a fixed ground anchor. It does not include a wing or bridle system, making it ideal for isolating and analyzing the behavior of the tethers themselves.
Arguments
set::Settings
: Configuration parameters.
Keywords
axial_stiffness::Vector{Float64}
: Predefined axial stiffness [N] for each tether.axial_damping::Vector{Float64}
: Predefined axial damping [Ns] for each tether.
Returns
SystemStructure
: A newSystemStructure
representing the 4-tether test system.
SymbolicAWEModels.create_simple_ram_sys_struct
— Functioncreate_simple_ram_sys_struct(set::Settings; axial_stiffness, axial_damping)
Create a simplified SystemStructure
for a ram-air kite with direct tether connections.
This model represents a kite with a flexible wing (2 deformable groups) but simplifies the bridle by connecting the four main tethers directly to the wing attachment points, omitting the complex pulley system. Each tether is modeled as a single segment.
Arguments
set::Settings
: Configuration parameters.
Keywords
axial_stiffness::Vector{Float64}
: Predefined axial stiffness [N] for each tether.axial_damping::Vector{Float64}
: Predefined axial damping [Ns] for each tether.
Returns
SystemStructure
: A newSystemStructure
representing the simplified model.
State Accessor Functions (Getters)
Use these functions to retrieve state information and calculated values from a model instance.
SymbolicAWEModels.winch_force
— FunctionReturns the winch force [N] for each winch.
SymbolicAWEModels.unstretched_length
— FunctionReturns the unstretched tether length [m] for each winch.
SymbolicAWEModels.tether_length
— FunctionReturns the current stretched tether length [m] for each tether.
Utility and Helper Functions
General helper functions for package management and setup.
SymbolicAWEModels.copy_examples
— Functioncopy_examples()
Copy all example scripts to the folder "examples" (it will be created if it doesn't exist).
SymbolicAWEModels.copy_bin
— Functioncopy_bin()
Copy the scripts createsysimage and run_julia to the folder "bin" (it will be created if it doesn't exist).