High-Level Functions
FLORIS
FLORIDyn.calcCt
— FunctioncalcCt(a::Number, _) -> Float64
calcCt(a::AbstractVector, _) -> AbstractVector
Calculate the thrust coefficient ct = 4a(1-a).
FLORIDyn.centerline!
— Functioncenterline!(deflection::AbstractMatrix,
states_op, states_t, states_wf, floris::Floris, d_rotor)
Compute the cross-wind wake deflection at the observation points in-place using the Gaussian wake model.
This function calculates the lateral (y) and vertical (z) deflection of the wake centerline due to yaw misalignment and other effects. The results are written directly into the provided deflection matrix without allocating temporary arrays.
Output Arguments
deflection::AbstractMatrix
(size n×2): Wake deflection components filled in-place- Column 1: Lateral deflection Δy
[m]
- Column 2: Vertical deflection Δz
[m]
(always zero in current implementation)
- Column 1: Lateral deflection Δy
Input Arguments
states_op::AbstractMatrix
(size n×k): Observation point states where n is number of points- Column 4 contains downstream distance in wake-aligned coordinates
[m]
- Column 4 contains downstream distance in wake-aligned coordinates
states_t::AbstractMatrix
: Turbine state matrix containing:- Column 1: Axial induction factor a
[-]
- Column 2: Yaw angle
[degrees]
- Column 3: Local turbulence intensity TI
[-]
- Column 1: Axial induction factor a
states_wf::AbstractMatrix
: Wind field state matrix containing:- Column 3: Ambient turbulence intensity TI₀
[-]
- Column 3: Ambient turbulence intensity TI₀
floris::Floris
: FLORIS model parameters for wake calculations (seeFloris
)d_rotor::Real
: Rotor diameter D[m]
Notes
- Only
states_op[:, 4]
(downstream distance) is used from the observation points - The function internally converts yaw angles from degrees to radians with sign correction
- Thrust coefficient is calculated from axial induction factor using
calcCt
See also: getVars!
FLORIDyn.discretizeRotor
— FunctiondiscretizeRotor(n_rp::Int) -> Tuple{Matrix{Float64}, Vector{Float64}}
Discretizes the rotor into n_rp
segments using the isocell algorithm.
Memoization: results are cached per thread. Repeated calls with the same n_rp
on the same thread reuse the cached arrays (no lock needed). Do not mutate the returned arrays, as they are shared within the thread.
Arguments
n_rp::Int
: The number of radial points to discretize the rotor into.
Returns
(m_rp, w)
where:m_rp::Matrix{Float64}
: Size(nC, 3)
; first column zeros, columns 2–3 are normalized coordinates in[-0.5, 0.5]
.w::Vector{Float64}
: Weights per cell that sum to approximately 1.
Notes
- Per-thread cache avoids contention; different threads may compute and hold their own cached copies for the same
n_rp
. - The isocell algorithm may not yield exactly
n_rp
cells but aims for a similar number. - For details, see: Masset et al. (2009) https://orbi.uliege.be/bitstream/2268/91953/1/massetisocellorbi.pdf
- The choice
N1 = 3
is used here; values of 4 or 5 are also viable.
FLORIDyn.init_states
— Functioninit_states(set::Settings, wf::WindFarm, wind::Wind, init_turb,
floris::Floris, sim::Sim) -> Tuple{Matrix, Matrix, Matrix}
Initialize the state arrays for wind farm simulation using the Gaussian wake model.
This function sets up the initial conditions for turbines, observation points, and wind field states based on the provided configuration parameters. It computes initial positions, wind conditions, and wake properties for each turbine in the wind farm.
Arguments
set::Settings
: Simulation settings containing configuration options for velocity, direction, and turbulence modelswf::WindFarm
: Wind farm object containing turbine positions, dimensions, and state arrays (seeWindFarm
)wind::Wind
: Wind conditions including velocity, direction, and turbulence intensity data (seeWind
)init_turb
: Initial turbine state parameters (axial induction factor, yaw angle, turbulence intensity)floris::Floris
: FLORIS model parameters for wake calculations (seeFloris
)sim::Sim
: Simulation parameters including time step and start time (seeSim
)
Returns
A tuple (states_op, states_t, states_wf)
containing:
states_op::Matrix
: Observation point states with 3D coordinates and wake positionsstates_t::Matrix
: Turbine states including control parameters and operational conditionsstates_wf::Matrix
: Wind field states with velocity, direction, and turbulence data
Description
The function performs the following initialization steps for each turbine:
- Retrieves wind field data (velocity, direction, turbulence intensity) based on the specified input methods
- Initializes wind field states at all observation points for the turbine
- Calculates downwind distances for wake coordinate system
- Sets initial turbine states from provided parameters
- Computes crosswind wake deflections using the centerline function
- Transforms coordinates from wake-relative to world coordinate system
- Updates observation point positions including turbine base and nacelle offsets
Notes
- Supports multiple wind input methods including interpolation, constant values, and random walk models
- Handles both 3D and 4D wind field configurations (with optional orientation data)
- Uses SOWFA coordinate system conventions for angle transformations
FLORIDyn.getVars!
— Function getVars!(sig_y, sig_z, x_0, delta, pc_y, pc_z, rps, c_t, yaw, ti, ti0, floris::Floris, d_rotor)
Compute Gaussian wake widths, deflection, potential-core radii, and onset distance at observation points, in-place.
Output Arguments
sig_y::AbstractVector{<:Real}
(length n): Lateral Gaussian width σ_y at each point [m]sig_z::AbstractVector{<:Real}
(length n): Vertical Gaussian width σ_z at each point [m]x_0::AbstractVector{<:Real}
(length n): Onset distance of the far-wake x₀ [m]delta::AbstractMatrix{<:Real}
(length n×2): Deflection components[Δy, Δz]
[m]pc_y::AbstractVector{<:Real}
(length n): Potential-core radius in y at each point [m]pc_z::AbstractVector{<:Real}
(length n): Potential-core radius in z at each point [m]
Input Arguments
rps::AbstractMatrix
(n×3): Observation points in wake-aligned frame; columns are[x_downstream, y_cross, z_cross]
[m]c_t::Union{Number,AbstractVector}
: Thrust coefficient Ct (scalar or length n) [-]yaw::Union{Number,AbstractVector}
: Yaw misalignment (scalar or length n) [rad]ti::Union{Number,AbstractVector}
: Local turbulence intensity TI at turbine (scalar or length n) [-]ti0::Union{Number,AbstractVector}
: Ambient turbulence intensity TI₀ (scalar or length n) [-]floris::Floris
: FLORIS Gaussian model parameters; seeFloris
d_rotor::Real
: Rotor diameter D [m]
Behavior
- Supports scalar or per-point values for
c_t
,yaw
,ti
,ti0
; scalars are broadcast to all points. - Uses
floris.k_a
,floris.k_b
,floris.alpha
,floris.beta
to compute per-pointx₀
,σ_y
,σ_z
, deflectionΔy
(hereΔz
is set to 0), and potential-core radiipc_y
,pc_z
. - No heap allocations beyond the provided outputs; results are written in-place and the function returns
nothing
.
Notes
- Only
rps[:, 1]
(downstream distance) is used by this implementation;rps[:, 2:3]
are ignored. delta
must have at least 2 columns; only columns 1:2 are written.- Units: distances in meters, angles in radians, intensities and
Ct
are dimensionless.
Example
n = size(RPs, 1)
sig_y = similar(RPs[:, 1])
sig_z = similar(RPs[:, 1])
x0 = similar(RPs[:, 1])
delta = zeros(n, 2)
pc_y = similar(RPs[:, 1])
pc_z = similar(RPs[:, 1])
getVars!(sig_y, sig_z, x0, delta, pc_y, pc_z, RPs, Ct, yaw, TI, TI0, floris, D)
Returns
nothing
— all results are written into the provided arrays.
FLORIDyn.getPower
— FunctiongetPower(wf::WindFarm, m::AbstractMatrix, floris::Floris, con::Con)
Calculate the power output of wind turbines in a wind farm simulation.
This function computes the power generated by wind turbines based on their operational states, wind conditions, and control settings. It accounts for yaw angle effects and optional yaw range constraints using hyperbolic tangent functions for smooth operational limits.
Arguments
wf::WindFarm
: Wind farm object containing turbine states, dimensions, and operational data with current axial induction factors and yaw angles (seeWindFarm
)m::Matrix
: Measurement or simulation data matrix where column 3 contains effective wind speeds at turbine locations [m/s]floris::Floris
: FLORIS model parameters containing air density, drivetrain efficiency, and power curve parameters (seeFloris
)con::Con
: Controller configuration object with yaw control settings and operational constraints (seeCon
)
Returns
P::Vector{Float64}
: Power output for each turbine in the wind farm [W]
Description
The function calculates power using the standard wind turbine power equation with yaw corrections:
P = 0.5 × ρ × A × Cp × U³ × η × cos(γ)^p_p × f_yaw_constraints
Where:
ρ
is air density [floris.airDen
] [kg/m³]A
is rotor swept areaπ × (D/2)²
[m²]Cp
is power coefficient calculated as4a(1-a)²
[-]U
is effective wind speed from column 3 of matrixm
[m/s]η
is drivetrain efficiency [floris.eta
] [-]γ
is yaw angle [rad]p_p
is yaw power exponent [floris.p_p
] [-]f_yaw_constraints
is optional yaw range constraint factor [-]
The yaw constraint factor is applied when con.tanh_yaw
is true:
f_yaw_constraints = [0.5 × tanh((γ_max - γ) × 50) + 0.5] ×
[-0.5 × tanh((γ_min - γ) × 50) + 0.5]
Notes
- Power coefficient is calculated from axial induction factor:
Cp = 4a(1-a)²
- Yaw effects reduce power output according to
cos(γ)^p_p
wherep_p
is typically 1.88 - Optional yaw range constraints use hyperbolic tangent functions with slope factor 50 for smooth transitions
- When
con.tanh_yaw
is enabled, power is smoothly constrained within [con.yawRangeMin
,con.yawRangeMax
] - The constraint functions approach step functions but provide smooth gradients for optimization
- Axial induction factors are extracted from
wf.States_T[wf.StartI, 1]
for current time step - Yaw angles are converted from degrees to radians internally
FLORIDyn.runFLORIS!
— FunctionrunFLORIS!(buffers::FLORISBuffers, set::Settings, location_t, states_wf, states_t, d_rotor,
floris, windshear::Union{Matrix, WindShear})
Execute the FLORIS (FLOw Redirection and Induction in Steady State) wake model simulation for wind farm analysis.
This is the main orchestrating function that coordinates the FLORIS wake model execution through a series of specialized sub-functions. It performs wake analysis using the Gaussian wake model to calculate velocity reductions, turbulence intensity additions, and effective wind speeds at turbine locations, accounting for wake interactions, rotor discretization, wind shear effects, and turbulence propagation.
Arguments
buffers::FLORISBuffers
: Pre-allocated buffer arrays for computation and output storage (seeFLORISBuffers
)set::Settings
: Simulation settings containing configuration options for wind shear modelinglocation_t
: Matrix of turbine positions [x, y, z] coordinates for each turbine [m]states_wf
: Wind field state matrix containing velocity, direction, and turbulence datastates_t
: Turbine state matrix with axial induction factors, yaw angles, and turbulence intensitiesd_rotor
: Vector of rotor diameters for each turbine [m]floris
: FLORIS model parameters containing wake model coefficients and rotor discretization settings (seeFloris
)windshear
: Wind shear profile data for vertical wind speed variation modeling, either a matrix or of typeWindShear
Returns
nothing
. Results are written into fields ofbuffers
:buffers.T_red_arr
: Velocity reduction factors for each turbinebuffers.T_aTI_arr
: Added turbulence intensity valuesbuffers.T_Ueff
: Effective wind speedsbuffers.T_weight
: Wake weighting factors
Implementation Structure
The function is implemented as a high-level orchestrator that calls specialized sub-functions:
prepare_rotor_points!
: Handles rotor discretization, scaling, rotation, and translationhandle_single_turbine!
: Optimized path for single turbine simulations (wind shear only)setup_computation_buffers!
: Initializes and configures computation buffers and viewscompute_wake_effects!
: Core wake computation loop for each upstream turbinecompute_final_wind_shear!
: Final wind shear and effective wind speed calculation
Computational Process
Single Turbine Case
For single turbine simulations, the function uses an optimized path that only calculates wind shear effects without wake interactions, providing significant performance benefits.
Multi-Turbine Wake Analysis
For multi-turbine wind farms, the function performs:
Rotor Point Preparation
- Discretizes rotor planes using the isocell algorithm
- Applies yaw rotation transformations and coordinate translations
- Handles both active turbines (d_rotor > 0) and placeholder turbines
Wake Interaction Calculations
For each upstream turbine affecting downstream turbines:
- Coordinate Transformations: Aligns coordinates with wind direction and turbine yaw
- Wake Variable Calculations: Computes wake expansion, deflection, and potential core dimensions
- Velocity Deficit Modeling: Applies Gaussian wake theory with yaw corrections
- Turbulence Enhancement: Calculates added turbulence using empirical correlations
- Superposition: Combines effects from multiple upstream turbines
Final Integration
- Applies vertical wind shear corrections to the downstream turbine
- Combines all wake effects with wind shear for final effective wind speed
Mathematical Models
The function implements state-of-the-art wake modeling based on:
- Gaussian Wake Theory: For velocity deficit calculations with yaw corrections
- Analytical Deflection Models: For wake steering in yawed conditions
- Empirical Turbulence Models: For wake-added turbulence intensity
- Linear Superposition: For combining multiple wake interactions
Performance Characteristics
- Computational Complexity: O(N²) for N turbines due to wake interactions
- Memory Efficiency: Uses pre-allocated buffers to avoid runtime allocations
- Optimization: Single turbine case bypasses multi-turbine calculations
- Vectorization: Leverages SIMD operations where possible
Notes
- Uses SOWFA (Simulator for Offshore Wind Farm Applications) coordinate conventions
- Supports both research and engineering applications for wind farm optimization
- Requires proper initialization of turbine states and wind field conditions
- Buffer sizes must be compatible with rotor discretization settings
References
- Bastankhah, M. and Porté-Agel, F. (2016). Experimental and theoretical study of wind turbine wakes in yawed conditions
- Niayifar, A. and Porté-Agel, F. (2016). Analytical modeling of wind farms: A new approach for power prediction
See Also
prepare_rotor_points!
: Rotor point preparation and transformationhandle_single_turbine!
: Single turbine optimization pathsetup_computation_buffers!
: Buffer initialization and setupcompute_wake_effects!
: Core wake interaction calculationscompute_final_wind_shear!
: Final wind shear integrationFLORISBuffers
: Buffer structure documentationFloris
: FLORIS model parameters
FLORIDyn.prepare_rotor_points!
— Functionprepare_rotor_points!(buffers::FLORISBuffers, location_t, states_t, d_rotor, floris::Floris)
Prepare rotor discretization points with scaling, rotation, and translation for the last turbine.
This function handles the rotor point discretization, scales them by the rotor diameter, applies yaw rotation, and translates them to the turbine location. The results are stored in the buffers to avoid allocations.
Arguments
buffers::FLORISBuffers
: Pre-allocated computation bufferslocation_t
: Turbine locations matrixstates_t
: Turbine states matrix (includes yaw angles)d_rotor
: Rotor diameter arrayfloris::Floris
: FLORIS model parameters
Returns
RPl
: View of transformed rotor pointsRPw
: Rotor point weights
Note
This function is private and intended for internal use only.
FLORIDyn.handle_single_turbine!
— Functionhandle_single_turbine!(buffers::FLORISBuffers, RPl, RPw, location_t, set::Settings, windshear, d_rotor)
Handle the special case when there is only one turbine in the simulation.
This function computes the wind shear reduction for a single turbine case and populates the appropriate buffer arrays. It returns early to avoid the multi-turbine wake calculations.
Arguments
buffers::FLORISBuffers
: Pre-allocated computation buffersRPl
: Rotor discretization pointsRPw
: Rotor point weightslocation_t
: Turbine locations matrixset::Settings
: Simulation settingswindshear
: Wind shear datad_rotor
: Rotor diameter array
Returns
nothing
(results stored in buffers)
Note
This function is private and intended for internal use only.
FLORIDyn.setup_computation_buffers!
— Functionsetup_computation_buffers!(buffers::FLORISBuffers, nRP::Int, nT::Int)
Initialize and setup computation buffers for multi-turbine wake calculations.
This function resizes output arrays and creates views of pre-allocated buffers to match the current rotor discretization size. It ensures all buffers are properly sized before the main computation loop.
Arguments
buffers::FLORISBuffers
: Pre-allocated computation buffersnRP::Int
: Number of rotor pointsnT::Int
: Number of turbines
Returns
- Tuple of buffer views for use in wake calculations
Note
This function is private and intended for internal use only.
FLORIDyn.compute_wake_effects!
— Functioncompute_wake_effects!(buffers::FLORISBuffers, views, iT::Int, RPl, RPw, location_t,
states_wf, states_t, d_rotor, floris::Floris, nRP::Int)
Compute wake effects for a single upstream turbine on the downstream turbine.
This function calculates the velocity deficit and added turbulence caused by turbine iT on the last turbine in the array. It handles coordinate transformations, wake variable calculations, and Gaussian wake modeling.
Arguments
buffers::FLORISBuffers
: Pre-allocated computation buffersviews
: Tuple of buffer views from setupcomputationbuffers!iT::Int
: Index of the upstream turbineRPl
: Rotor discretization pointsRPw
: Rotor point weightslocation_t
: Turbine locations matrixstates_wf
: Wind farm states matrixstates_t
: Turbine states matrixd_rotor
: Rotor diameter arrayfloris::Floris
: FLORIS model parametersnRP::Int
: Number of rotor points
Returns
nothing
. Results are written into fields ofbuffers
:buffers.T_red_arr
: Velocity reduction factors for each turbinebuffers.T_aTI_arr
: Added turbulence intensity valuesbuffers.T_Ueff
: Effective wind speedsbuffers.T_weight
: Wake weighting factors
Note
This function is private and intended for internal use only.
FLORIDyn.compute_final_wind_shear!
— Functioncompute_final_wind_shear!(buffers::FLORISBuffers, RPl, RPw, location_t, set::Settings,
windshear, tmp_RPs_r, states_wf)
Compute final wind shear effects and effective wind speed for the last turbine.
This function calculates the wind shear reduction for the downstream turbine and computes the final effective wind speed by combining all wake effects and wind shear.
Arguments
buffers::FLORISBuffers
: Pre-allocated computation buffersRPl
: Rotor discretization pointsRPw
: Rotor point weightslocation_t
: Turbine locations matrixset::Settings
: Simulation settingswindshear
: Wind shear datatmp_RPs_r
: Temporary buffer for rotor point calculationsstates_wf
: Wind farm states matrix
Returns
nothing
. Results are written into fields ofbuffers
:buffers.T_red_arr
: Velocity reduction factors for each turbinebuffers.T_aTI_arr
: Added turbulence intensity valuesbuffers.T_Ueff
: Effective wind speedsbuffers.T_weight
: Wake weighting factors
Note
This function is private and intended for internal use only.
FLORIDyn
FLORIDyn.initSimulation
— FunctioninitSimulation(wf::Union{Nothing, WindFarm}, sim::Sim) -> Union{Nothing, WindFarm}
Initialize or load a wind farm simulation state based on simulation settings.
This function handles the initialization phase of a wind farm simulation by either saving the current initialized state to disk or loading a previously saved state, depending on the simulation configuration.
Arguments
wf::Union{Nothing, WindFarm}
: Wind farm object containing the initialized simulation state, orNothing
if no state is available. SeeWindFarm
sim::Sim
: Simulation configuration object containing initialization settings and file paths. SeeSim
Returns
wf::Union{Nothing, WindFarm}
: The wind farm state, either the original input state (for "init" mode) or a loaded state from disk (for "load" mode)
Behavior
The function operates in two modes based on sim.init
:
"init" Mode
- Uses the provided wind farm state as-is
- If
sim.save_init_state
istrue
, saves the current state to"T_init.jld2"
in the specified data directory - Logs the save operation for user feedback
"load" Mode
- Attempts to load a previously saved wind farm state from
"T_init.jld2"
- Falls back to the provided state if loading fails (with warning)
- Handles file I/O errors gracefully
File Operations
- Save path:
$(sim.path_to_data)/T_init.jld2
- Format: JLD2 binary format for efficient Julia object serialization
- Error handling: Loading failures produce warnings but do not halt execution
Notes
- The function is case-insensitive for the initialization mode string
- File operations use the path specified in
sim.path_to_data
- Loading errors are caught and logged as warnings, allowing simulation to proceed with the original state
- This mechanism enables reproducible simulations by preserving and reusing initial conditions
FLORIDyn.findTurbineGroups
— FunctionfindTurbineGroups(wf::WindFarm, floridyn::FloriDyn) -> Vector{Vector{Int64}}
Determine wake interaction dependencies between turbines in a wind farm.
This function analyzes the spatial relationships between turbines to identify which turbines are affected by the wakes of upstream turbines. It uses coordinate transformations to the wind-aligned reference frame and geometric criteria to determine wake interactions.
Arguments
wf::WindFarm
: Wind farm object containing turbine positions, observation points, and wind field states. SeeWindFarm
floridyn::FloriDyn
: FLORIDyn model parameters containing wake interaction thresholds. SeeFloriDyn
Returns
vv_dep::Vector{Vector{Int64}}
: A vector of vectors wherevv_dep[i]
contains the indices of all turbines that affect turbinei
through wake interactions. Each inner vector lists the upstream turbine indices that influence the wake conditions at the corresponding turbine.
Algorithm
- Coordinate Transformation: For each turbine pair, transforms coordinates to a wind-aligned frame using the closest observation point
- Wake Zone Detection: Applies geometric criteria to determine if a downstream turbine lies within the wake zone:
- Upstream extent:
r₁[1] ≥ -uw × D[iT]
(allowing for slight upstream influence) - Downstream extent:
r₁[1] ≤ dw × D[iT]
(wake extends downstream) - Lateral extent:
|r₁[2]| ≤ cw × D[iT]
(wake width constraint)
- Upstream extent:
- Dependency Matrix: Constructs a boolean dependency matrix and extracts indices for each turbine
Mathematical Description
The wake interaction criteria are evaluated in the wind-aligned coordinate system:
r₁ = R(φ) × (rₒₚ - rₜᵤᵣᵦ)
where:
R(φ)
is the rotation matrix for wind direction angleφ
rₒₚ
is the position of the closest observation point from the upstream turbinerₜᵤᵣᵦ
is the position of the downstream turbine being evaluated
Notes
- The function uses the closest observation point from each upstream turbine to determine wind direction
- Wake zones are defined as multiples of rotor diameter using the FLORIDyn parameters
- Self-interaction (turbine affecting itself) is explicitly excluded
- The coordinate transformation accounts for the SOWFA wind direction convention
FLORIDyn.prepareSimulation
— FunctionprepareSimulation(set::Settings, wind::Wind, con::Con, floridyn::FloriDyn,
floris::Floris, turbProp, sim::Sim) -> (WindFarm, Wind, Sim, Con, Floris)
Prepares the simulation environment for a wind farm analysis using the provided settings and parameters.
Arguments
set::Settings
: Simulation settings containing configuration options.wind::Wind
: Wind conditions or wind field data. See:Wind
con::Con
: Controller parameters of the turbines. See:Con
floridyn::FloriDyn
: Parameters specific to the FLORIDyn model. See:FloriDyn
floris::Floris
: Parameters specific to the FLORIS model. See:Floris
turbProp
: Properties of the turbines involved in the simulation.sim::Sim
: Simulation-specific parameters or state. See:Sim
Arguments that get modified
wind
: Updated with wind velocity, direction, turbulence intensity, and shear profile.con
: Updated with yaw data.sim
: Updated with the number of simulation steps.floris
: May include additional parameters for the FLORIS model.
Returns
- Returns the tuple
(wf, wind, sim, con, floris)
where:wf
: Wind farm struct containing turbine states and positions. See:WindFarm
wind
: Updated wind conditions.sim
: Updated simulation parameters.con
: Updated controller parameters.floris
: Parameters for the FLORIS model.
FLORIDyn.perturbationOfTheWF!
— FunctionperturbationOfTheWF!(wf::WindFarm, wind::Wind) -> Nothing
Apply stochastic perturbations to the wind field states in-place.
This function adds Gaussian noise to the wind field parameters to model measurement uncertainty or natural variability in wind conditions. The perturbations are applied conditionally based on the wind perturbation configuration and are added directly to the wind farm state matrix.
Input/ Output Arguments
wf::WindFarm
: Wind farm struct containing the state matrixStates_WF
to be perturbed
Input Arguments
wind::Wind
: Wind configuration struct containing perturbation settings. SeeWind
Returns
nothing
: The function modifies the wind farm state in-place
Behavior
The function applies independent Gaussian perturbations to three wind field parameters:
Velocity Perturbation
- Condition:
wind.perturbation.vel == true
- Target: Column 1 of
wf.States_WF
(wind velocity [m/s]) - Noise:
wind.perturbation.vel_sigma * randn(nOP × nT)
Direction Perturbation
- Condition:
wind.perturbation.dir == true
- Target: Column 2 of
wf.States_WF
(wind direction [degrees]) - Noise:
wind.perturbation.dir_sigma * randn(nOP × nT)
Turbulence Intensity Perturbation
- Condition:
wind.perturbation.ti == true
- Target: Column 3 of
wf.States_WF
(turbulence intensity [-]) - Noise:
wind.perturbation.ti_sigma * randn(nOP × nT)
Mathematical Description
For each enabled perturbation type, the function applies:
States_WF[:, col] += σ × N(0,1)
where:
σ
is the standard deviation for the specific parameterN(0,1)
is standard normal random noise with dimensions(nOP × nT)
nOP
is the number of observation points per turbinenT
is the total number of turbines
Notes
- The function uses in-place modification (indicated by the
!
suffix) - Perturbations are applied independently to each observation point and turbine
- The random noise follows a standard normal distribution scaled by the respective sigma values
- Only enabled perturbation types (based on boolean flags) are applied
FLORIDyn.setUpTmpWFAndRun!
— FunctionsetUpTmpWFAndRun!(ub::UnifiedBuffers, wf::WindFarm, set::Settings, floris::Floris,
wind::Wind) -> Nothing
Non-allocating version that uses a unified buffer struct for wind farm calculations.
This function performs wind farm wake calculations while avoiding memory allocations by reusing pre-allocated buffer arrays from a UnifiedBuffers
struct. This is particularly important for parallel execution and performance-critical loops where garbage collection overhead needs to be minimized.
Buffer Arguments
ub::UnifiedBuffers
: Unified buffer struct containing all pre-allocated arraysub.M_buffer
: Pre-allocated buffer for results matrix (size: nT × 3)ub.iTWFState_buffer
: Buffer for turbine wind field stateub.tmp_Tpos_buffer
: Buffer for temporary turbine positionsub.tmp_WF_buffer
: Buffer for temporary wind field statesub.tmp_Tst_buffer
: Buffer for temporary turbine statesub.dists_buffer
: Buffer for distance calculationsub.plot_WF_buffer
: Buffer for plotting wind field dataub.plot_OP_buffer
: Buffer for plotting operating point data
Output Arguments
ub.M_buffer
: Pre-allocated buffer for results matrix (size: nT × 3)wf.Weight
: Sets wake weight factors for each turbine from FLORIS calculationswf.red_arr
: Updates wake reduction factors between turbines (wake interference matrix)
Input/ Output Arguments
wf::WindFarm
: Wind farm object containing turbine data
Input Arguments
set::Settings
: Settings object containing simulation parametersfloris::Floris
: FLORIS model parameters for wake calculationswind::Wind
: Wind field configuration
Returns
- nothing: The function modifies
ub.M_buffer
,wf.Weight
, andwf.red_arr
in-place, storing the results of the wind farm calculations
Performance Notes
- Uses in-place operations to minimize memory allocations
- Buffers must be pre-sized correctly for the specific wind farm configuration
- Thread-safe when each thread uses its own set of buffers
FLORIDyn.interpolateOPs!
— FunctioninterpolateOPs!(unified_buffers::UnifiedBuffers, intOPs::Vector{Matrix{Float64}},
wf::WindFarm) -> Nothing
Compute interpolation weights and indices for observation points affecting each turbine using a unified buffer.
This function performs interpolation calculations while avoiding memory allocations by reusing pre-allocated buffer arrays from a unified buffer struct. This is critical for performance when called repeatedly in loops, such as in flow field calculations.
Buffer Arguments
unified_buffers::UnifiedBuffers
: Unified buffer struct containing pre-allocated arrays including:dist_buffer
: Buffer for distance calculations (length ≥ wf.nOP)sorted_indices_buffer
: Buffer for sorting indices (length ≥ wf.nOP)
Output Arguments
intOPs::Vector{Matrix{Float64}}
: Pre-allocated vector of matrices to store interpolation results
Input Arguments
wf::WindFarm
: Wind farm object containing turbine positions and observation point data
Returns
- nothing: The function modifies
intOPs
in-place, storing the interpolation results for each turbine
Example
# Create unified buffers
unified_buffers = create_unified_buffers(wf)
# Pre-allocate interpolation matrices
intOPs = [zeros(length(wf.dep[iT]), 4) for iT in 1:wf.nT]
# Non-allocating interpolation
interpolateOPs!(unified_buffers, intOPs, wf)
FLORIDyn.iterateOPs!
— FunctioniterateOPs!(iterate_mode::IterateOPs_model, wf::WindFarm, sim::Sim, floris::Floris,
floridyn::FloriDyn, buffers::IterateOPsBuffers) -> Nothing
Advance observation points through the wind field using the specified iteration strategy.
This function family implements different algorithms for moving observation points (OPs) through space and time, which is essential for accurate wake propagation modeling in wind farm simulations. The choice of iteration method affects computational efficiency, numerical stability, and physical accuracy.
Summary
The function modifies the following WindFarm fields:
wf.States_OP
: Updates observation point positions and states through temporal advancementwf.States_T
: Updates turbine states through circular shifting and temporal evolutionwf.States_WF
: Updates wind field states through circular shifting and temporal evolution
Input/ Output Arguments
wf::WindFarm
: Wind farm object containing turbine and observation point data
Input Arguments
iterate_mode::IterateOPs_model
: Iteration strategy (e.g.,IterateOPs_basic
,IterateOPs_average
)sim::Sim
: Simulation configuration with time-stepping parametersfloris::Floris
: FLORIS model parameters for wake calculationsfloridyn::FloriDyn
: FLORIDyn model parameters for wake dynamicsbuffers::IterateOPsBuffers
: Pre-allocated buffers for allocation-free execution
Algorithm Overview
- State Preservation: Save initial turbine observation point states
- Downwind Advection: Move OPs downstream based on local wind velocity
- Crosswind Deflection: Apply wake-induced lateral deflection using centerline calculations
- Coordinate Transformation: Convert to world coordinates using wind direction
- Temporal Advancement: Perform circular shifting to advance time steps
- Spatial Reordering: Maintain downstream position ordering of observation points
Available Methods
iterateOPs!(::IterateOPs_basic, ...)
: Basic time-stepping with simple advection
Notes
- Different iteration strategies provide trade-offs between accuracy and computational cost
FLORIDyn.angSOWFA2world
— FunctionangSOWFA2world(deg_SOWFA) -> Float64
Convert wind direction angle from SOWFA convention to world coordinate system.
This function performs coordinate transformation between different wind direction conventions used in wind farm simulations. SOWFA (Simulator fOr Wind Farm Applications) uses a different angular reference system than the standard world coordinate system used in calculations.
Arguments
deg_SOWFA::Real
: Wind direction angle in SOWFA convention [degrees]
Returns
rad_World
: Wind direction angle in world coordinate system [radians]
Coordinate System Conversion
The transformation follows the relationship:
θ_world = 270° - θ_SOWFA
SOWFA Convention
- Wind direction angles are defined clockwise from a reference direction
World Convention
- Wind direction angles are defined counterclockwise for mathematical calculations
- Standard convention used in wake models and analytical computations
Mathematical Description
The conversion process:
- Angular transformation:
deg_World = 270 - deg_SOWFA
- Unit conversion:
rad_World = deg2rad(deg_World)
The 270° offset accounts for the difference between clockwise (SOWFA) and counterclockwise (world) angular conventions.
Examples
# Convert 90° SOWFA direction to world coordinates
world_angle = angSOWFA2world(90.0) # Returns 3.141592... (180° in radians)
# Convert 0° SOWFA direction
world_angle = angSOWFA2world(0.0) # Returns 4.712388... (270° in radians)
Notes
- The function handles the sign convention difference between coordinate systems
- Output is always in radians for use in trigonometric calculations
- This transformation is essential for proper wake modeling in wind farm simulations
- The 270° offset ensures proper alignment between SOWFA and mathematical conventions
FLORIDyn.runFLORIDyn
— FunctionrunFLORIDyn(plt, set::Settings, wf::WindFarm, wind::Wind, sim, con, vis, floridyn, floris;
rmt_plot_fn=nothing, msr=VelReduction) -> (WindFarm, DataFrame, Matrix)
Main entry point for the FLORIDyn closed-loop simulation.
Arguments
plt
: Plot object for live visualization during simulationset::Settings
: Simulation settings and configuration parameters.wf::WindFarm
: See: WindFarm simulation state, including turbine and wind farm states.wind::Wind
: See: Wind field settings.sim::Sim
: Simulation state or configuration object. See:Sim
con::Con
: Controller object or control parameters. See:Con
vis::Vis
: Visualization settings controlling online plotting and animation. See:Vis
floridyn::FloriDyn
: Parameters specific to the FLORIDyn model. See:FloriDyn
floris::Floris
: Parameters specific to the FLORIS model. See:Floris
Keyword Arguments
rmt_plot_fn
: Optional remote plotting function for intermediate simulation results. When provided, this function is called remotely (using@spawnat 2
) to plot flow field visualization on a separate worker process. The function should accept parameters(wf, X, Y, Z, vis, t_rel; msr=VelReduction)
wherewf
is the wind farm state,X
,Y
,Z
are flow field coordinates and velocities,vis
contains visualization settings, andt_rel
is the relative simulation time. Defaults tonothing
for local plotting.msr
: Measurement type for velocity reduction calculations. Defaults toVelReduction
.
Returns
A tuple (wf, md, mi)
containing:
wf::WindFarm
: Updated simulation state with final turbine positions, wind field states, and observation point datamd::DataFrame
: Measurement data with columns::Time
: Simulation time steps:ForeignReduction
: Wind speed reduction factors (%) due to wake effects from other turbines:AddedTurbulence
: Additional turbulence intensity (%) induced by upstream turbines:EffWindSpeed
: Effective wind speed (m/s) at each turbine after wake effects:FreeWindSpeed
: Free-stream wind speed (m/s) without wake interference:PowerGen
: Generated electrical power (MW) for each turbine
mi::Matrix
: Interaction matrix combining time data with turbine-to-turbine wake interaction coefficients for each simulation step
Description
Runs a closed-loop wind farm simulation using the FLORIDyn and FLORIS models, applying control strategies and updating turbine states over time.
FLORIDyn.run_floridyn
— Functionrun_floridyn(plt, set, wf, wind, sim, con, vis,
floridyn, floris; msr=VelReduction) -> (WindFarm, DataFrame, Matrix)
Unified function that automatically handles both multi-threading and single-threading modes for running FLORIDyn simulations with appropriate plotting callbacks.
Arguments
plt
: PyPlot instance, usually provided by ControlPlotsset
: Settings object. See: Settingswf
: WindFarm struct. These are work arrays, not persistent objects. See: WindFarmwind
: Wind field input settings. See: Windsim
: Simulation settings. See: Simcon
: Controller settings. See: Convis
: Visualization settings. See: Visfloridyn
: FLORIDyn model struct. See: FloriDynfloris
: Floris model struct. See: Florismsr
: Measurement index for online flow field plotting (VelReduction, AddedTurbulence or EffWind). Default VelReduction. See: MSR
Returns
- Tuple (wf, md, mi): WindFarm, measurement data, and interaction matrix
Visualization
FLORIDyn.create_thread_buffers
— Functioncreate_thread_buffers(wf::WindFarm, nth::Int, floris::Floris) -> ThreadBuffers
Create thread-local buffers for parallel flow field computation with FLORIS parameters.
This function pre-allocates all necessary data structures for each thread to avoid race conditions and memory allocations during the parallel computation loop.
Arguments
wf::WindFarm
: Original wind farm object to use as templatenth::Int
: Number of threads to create buffers forfloris::Floris
: FLORIS parameters for creating proper FLORIS buffers
Returns
ThreadBuffers
: Struct containing all thread-local buffers
Performance Notes
- Each thread gets its own copy of the WindFarm structure
- Pre-allocates all arrays to minimize allocations during computation
- Sets up dependency structure for virtual turbines at grid points
create_thread_buffers(wf::WindFarm, nth::Int) -> ThreadBuffers
Create thread-local buffers for parallel flow field computation.
This function pre-allocates all necessary data structures for each thread to avoid race conditions and memory allocations during the parallel computation loop.
Arguments
wf::WindFarm
: Original wind farm object to use as templatenth::Int
: Number of threads to create buffers for
Returns
ThreadBuffers
: Struct containing all thread-local buffers
Performance Notes
- Each thread gets its own copy of the WindFarm structure
- Pre-allocates all arrays to minimize allocations during computation
- Sets up dependency structure for virtual turbines at grid points
FLORIDyn.getMeasurements
— FunctiongetMeasurements(buffers::ThreadBuffers, mx::Matrix, my::Matrix, nM::Int, zh::Real,
wf::WindFarm, set::Settings, floris::Floris, wind::Wind) -> Array{Float64,3}
Calculate flow field measurements at specified grid points by treating them as virtual turbines.
This function computes flow field properties (velocity reduction, added turbulence, effective wind speed) at grid points by creating virtual turbines at each location and running the FLORIS wake model. Each grid point is treated as a turbine that depends on all real turbines in the wind farm, allowing wake effects to be captured in the flow field visualization.
Arguments
buffers::ThreadBuffers
: Pre-allocated thread-local buffers created withcreate_thread_buffers
; for Julia 1.12 usecreate_thread_buffers(wf, nthreads() + 1, floris)
; for single-thread usecreate_thread_buffers(wf, 1, floris)
mx::Matrix
: X-coordinates of grid points (m)my::Matrix
: Y-coordinates of grid points (m)nM::Int
: Number of measurements to compute (typically 3)zh::Real
: Hub height for measurements (m)wf::WindFarm
: Wind farm object containing turbine data. See:WindFarm
wf.nT
: Number of real turbineswf.StartI
: Starting indices for turbine datawf.posBase
,wf.posNac
: Turbine positionswf.States_*
: Turbine state matrices
set::Settings
: Settings object containing simulation parameters. See:Settings
floris::Floris
: FLORIS model parameters for wake calculations. See:Floris
wind::Wind
: Wind field configuration. See:Wind
Returns
mz::Array{Float64,3}
: 3D array of measurements with dimensions(size(mx,1), size(mx,2), nM)
mz[:,:,1]
: Velocity reductionmz[:,:,2]
: Added turbulence intensitymz[:,:,3]
: Effective wind speed
Algorithm
For each grid point:
- Creates a temporary wind farm with all original turbines plus one virtual turbine at the grid point
- Sets the virtual turbine to depend on all real turbines (to capture wake effects)
- Runs the FLORIS simulation to compute wake-affected flow properties
- Extracts the result for the virtual turbine position
Performance Notes
- Multi-threaded implementation using
@threads
for parallel processing of grid points when more than one buffer is provided - With a single buffer (
length(buffers.thread_buffers) == 1
), runs in a single-threaded loop - Each grid point requires a full wind farm simulation, so computation time scales with grid size
- Uses thread-local buffers created by
create_thread_buffers
to avoid race conditions - On Julia 1.12 create
nthreads() + 1
buffers to accommodate thread indexing
Example
# Create a 10x10 grid from 0 to 1000m
x_range = 0:100:1000
y_range = 0:100:1000
mx = repeat(collect(x_range)', length(y_range), 1)
my = repeat(collect(y_range), 1, length(x_range))
# Calculate 3 measurements at 90m hub height (single-thread)
buffers = create_thread_buffers(wind_farm, 1, floris_model)
mz = getMeasurements(buffers, mx, my, 3, 90.0, wind_farm, settings, floris_model, wind_config)
# Extract effective wind speed field
wind_speed_field = mz[:, :, 3]
See Also
calcFlowField
: Higher-level function that uses this to create complete flow field datasetUpTmpWFAndRun!
: Underlying simulation function used for each grid point
FLORIDyn.calcFlowField
— FunctioncalcFlowField(set::Settings, wf::WindFarm, wind::Wind, floris::Floris; plt=nothing)
Generate full flow field plot data by calculating measurements across a grid.
This function creates a rectangular grid over the wind farm domain and calculates flow field properties at each grid point by treating them as virtual turbines. The computation can be performed in parallel if set.threading
is true.
Arguments
set::Settings
: Settings object containing simulation parametersset.threading
: If true, uses multi-threaded computation with@threads
set.parallel
: If true, enables parallel-specific optimizations
wf::WindFarm
: Wind farm object containing turbine datawind::Wind
: Wind field configurationfloris::Floris
: FLORIS model parameters
Keyword Arguments
plt=nothing
: Plot object for garbage collection control. If provided andset.parallel
is true, automatically callsplt.GC.enable(false)
before multithreading andplt.GC.enable(true)
after completion to prevent PyCall-related segmentation faults during parallel execution with ControlPlots loaded. To take full advantage of multithreading, executed the plotting in a separate process.vis=nothing
: Visualization configuration object containing field limits and resolution settings. If provided, usesvis.field_limits_min
,vis.field_limits_max
, andvis.field_resolution
to define the computational grid. If not provided, defaults to domain [0,0,0] to [3000,3000,400] meters with 20m resolution.
Returns
Z::Array{Float64,3}
: 3D array of flow field measurements with dimensions(ny, nx, 3)
Z[:,:,1]
: Velocity reduction factorZ[:,:,2]
: Added turbulence intensityZ[:,:,3]
: Effective wind speed (m/s)
X::Matrix{Float64}
: X-coordinate grid (m)Y::Matrix{Float64}
: Y-coordinate grid (m)
Notes
- Grid resolution and domain are configurable via the
vis
parameter, or use default values for backward compatibility - Hub height is taken from the first turbine in the wind farm
Example
# Calculate flow field with threading and GC control
set.threading = true
Z, X, Y = calcFlowField(set, wf, wind, floris; plt)
# Extract velocity reduction field
velocity_reduction = Z[:, :, 1]
# Extract effective wind speed field
wind_speed = Z[:, :, 3]
See Also
getMeasurements
: Function used internally to compute the flow fieldplotFlowField
: Visualization function for the generated data
FLORIDyn.plotFlowField
— FunctionplotFlowField(state::Union{Nothing, PlotState}, plt, wf, mx, my, mz, vis::Vis, t=nothing;
msr::MSR=EffWind, fig=nothing)
Plot a 2D contour of the flow field data with support for animation.
Arguments
state::Union{Nothing, PlotState}
: Animation state object. Passnothing
for the first call, then pass the returned state for subsequent calls to maintain the same figure and layout.plt
: Plotting package (e.g., ControlPlots.plt)wf
: Wind farm object containing turbine datamx::Matrix
: X-coordinate gridmy::Matrix
: Y-coordinate gridmz::Array{Float64,3}
: 3D array of measurements with dimensions (rows, cols, nM)vis::Vis
: Visualization settings including save options and color scale parameterst
: Time value for display in the plot title or annotationsmsr::MSR
: Which measurement to plot. See: MSRvis.unit_test::Bool
: Whether to automatically close plots for testing.
Returns
state::PlotState
: Updated or newly created plot state for use in subsequent calls
Description
This function supports creating animations by maintaining plot state across multiple calls:
First Call (state = nothing)
- Creates new figure, axes, colorbar, and all plot elements
- Initializes and returns a PlotState object
Subsequent Calls (state = PlotState)
- Updates existing contour data, turbine positions, and observation points
- Reuses the same figure and layout for smooth animation
Animation Example
using ControlPlots
# Initialize state (first frame)
state = nothing
vis = Vis(online=true, save=true) # Enable saving to video folder
for t in time_steps
Z, X, Y = calcFlowField(settings, wind_farm, wind, floris)
state = plotFlowField(state, plt, wind_farm, X, Y, Z, vis, t; msr=EffWind)
plt.pause(0.01) # Small delay for animation
end
Notes
- The function automatically handles coordinate system transformations for turbine orientations
- Observation points are displayed as white scatter points for reference
- Color scales are kept consistent across animation frames when using the same measurement type
- The time parameter
t
can be used for title updates or time annotations - When
vis.save=true
, plots are saved as PNG files to thevideo/
directory - Saved filenames include measurement type and time information (e.g.,
velocity_reduction_t0120s.png
) - The
video/
directory is automatically created if it doesn't exist - This function requires a plotting package like ControlPlots.jl to be loaded and available as
plt
plotFlowField(plt, wf, mx, my, mz, vis, t=nothing; msr=EffWind, fig=nothing)
Compatibility method for the original plotFlowField interface.
This method provides backward compatibility by calling the new state-based version with state=nothing
, effectively creating a single plot without animation support.
Arguments
plt
: Plotting package (e.g., ControlPlots.plt)wf
: Wind farm object containing turbine datamx::Matrix
: X-coordinate gridmy::Matrix
: Y-coordinate gridmz::Array{Float64,3}
: 3D array of measurements with dimensions (rows, cols, nM)vis::Vis
: Visualization settings including save options and color scale parameterst
: Time value for display in the plot title or annotationsmsr::MSR
: Which measurement to plot. See: MSRvis.unit_test::Bool
: Whether to automatically close plots for testing.
Returns
nothing
: For compatibility with the original interface
Note
This method is provided for backward compatibility. For animation support, use the new interface with explicit state management.
FLORIDyn.plot_flow_field
— Functionplot_flow_field(wf, X, Y, Z, vis; msr=VelReduction, plt=nothing,
fig=nothing) -> Nothing
High-level plotting function that automatically dispatches to either parallel or sequential plotting based on the number of available threads and processes.
Arguments
wf
: WindFarm objectX
,Y
,Z
: Flow field coordinate arraysvis
: Visualization settingsmsr
: Measurement type, see: MSRplt
: Matplotlib PyPlot instance (only used in sequential mode)fig
: Figure name (optional)
Returns
- nothing
FLORIDyn.plotMeasurements
— FunctionplotMeasurements(plt, wf::WindFarm, md::DataFrame, vis::Vis;
separated=false, msr=VelReduction, pltctrl=nothing) -> Nothing
Plot measurement data from FLORIDyn simulation results.
Arguments
plt
: Plotting package (e.g., PyPlot, which is exported from ControlPlots); nothing for remote plottingwf::WindFarm
: Wind farm object with fieldnT
(number of turbines). SeeWindFarm
md::DataFrame
: Measurements DataFrame containing time series data with columns:Time
: Simulation time [s]ForeignReduction
: Foreign reduction [%] (for VelReduction)AddedTurbulence
: Added turbulence [%] (for AddedTurbulence)EffWindSpeed
: Effective wind speed [m/s] (for EffWind)
vis::Vis
: Visualization settings including unit_test parameter. SeeVis
separated::Bool
: Whether to use separated subplot layout (default: false)msr::MSR
: Measurement type to plot, see: MSR. Options: VelReduction, AddedTurbulence, EffWindpltctrl
: ControlPlots module instance (optional, used for large turbine count plotting); nothing for remote plotting
Returns
nothing
Description
This function creates time series plots of measurement data from FLORIDyn simulations. It handles:
- Time normalization by subtracting the start time
- Multiple measurement types (velocity reduction, added turbulence, effective wind speed)
- Automatic layout selection based on turbine count
- Both separated (subplot) and combined plotting modes
Plotting Modes
- Separated mode (
separated=true
): Creates individual subplots for each turbine- For ≤9 turbines: Traditional subplot grid
- For >9 turbines: Uses helper function with grouped subplots via
plot_x
- Combined mode (
separated=false
): Plots all turbines on a single figure with different colors
Examples
using ControlPlots # For single-threaded plotting
# Plot velocity reduction for all turbines in combined mode
plotMeasurements(plt, wind_farm, measurements_df, vis; msr=VelReduction)
# Plot added turbulence with separated subplots
plotMeasurements(plt, wind_farm, measurements_df, vis; separated=true, msr=AddedTurbulence)
# Plot effective wind speed with ControlPlots module for large farms
plotMeasurements(plt, wind_farm, measurements_df, vis; separated=true, msr=EffWind, pltctrl=ControlPlots)
See Also
plotFlowField
: For flow field visualizationgetMeasurements
: For generating measurement dataprepare_large_plot_inputs
: Helper for large turbine countsplot_x
: Multi-subplot plotting function
FLORIDyn.plot_measurements
— Functionplot_measurements(wf, md, vis; separated=true, msr=VelReduction, plt=nothing) -> Nothing
High-level measurements plotting function that automatically dispatches to either parallel or sequential plotting based on the number of available threads and processes.
Arguments
wf
: WindFarm objectmd
: Measurement datavis
: Visualization settingsseparated
: Whether to use separated subplotsmsr
: Measurement type, see: MSRplt
: Matplotlib PyPlot instance (only used in sequential mode)
Returns
- nothing
See Also
plotMeasurements
: The underlying plotting function used in sequential mode
FLORIDyn.plot_x
— Functionplot_x(times, plot_data...; ylabels=nothing, labels=nothing, fig="Wind Direction",
xlabel="rel_time [s]", ysize=10, bottom=0.02, legend_size=nothing, pltctrl=nothing,
loc=nothing) -> Nothing
High-level time series plotting function that automatically dispatches to either parallel or sequential plotting based on the number of available threads and processes.
Arguments
times
: Time vector for x-axisplot_data...
: Variable number of data arrays to plotylabels
: Labels for y-axes (optional)labels
: Labels for subplots (optional)fig
: Figure title (default: "Wind Direction")xlabel
: X-axis label (default: "rel_time [s]")ysize
: Size of the Y-axis labels in points (default: 10)bottom
: Bottom margin (default: 0.02)legend_size
: Legend font size in points (optional)pltctrl
: ControlPlots instance (only used in sequential mode)
Returns
- nothing
Description
When running with multiple threads and processes, it uses remote plotting capabilities via rmt_plotx
ONLY if no local pltctrl
is provided. If a pltctrl
argument is supplied (e.g. during unit tests with a mock), it forces the sequential path so tests can observe side-effects without needing remote worker setup.
Example
plot_x(times, data1, data2; ylabels=["Turbine 1", "Turbine 2"],
labels=["Wind Speed", "Power"], legend_size=8, pltctrl=pltctrl)
See Also
plotx
: The underlying plotting function from ControlPlots used in sequential mode
FLORIDyn.prepare_large_plot_inputs
— Functionprepare_large_plot_inputs(wf, md, data_column, ylabel) -> (times, plot_data, turbine_labels, subplot_labels)
Prepare grouped plotting inputs for large numbers of turbines when separated=true
.
Returns
times::Vector{Float64}
: Time vector (one entry per recorded time step)plot_data::Vector{Any}
: Each element is either a Vector (single line subplot) or Vector{Vector{Float64}} (multiple lines)turbine_labels::Vector{String}
: Y-axis labels per subplot (same ylabel repeated)subplot_labels::Vector{Vector{String}}
: Line labels per subplot
FLORIDyn.close_all
— Functionclose_all(plt)
Close all matplotlib figure windows.
This function automatically dispatches to either parallel or sequential plotting based on the number of available threads and processes.
Arguments
plt
: Matplotlib PyPlot instance (only used in sequential mode)
Description
When running with multiple threads and processes, it uses remote plotting capabilities to close all figures on the remote worker. Otherwise, it directly calls plt.close("all")
to close all figures in the current process.
Video Creation
FLORIDyn.cleanup_video_folder
— Functioncleanup_video_folder() -> Nothing
Clean up existing PNG files in the video folder before creating new videos.
Description
This function removes all PNG files from the "video" directory to ensure a clean slate before generating new video frames. It is typically called before running simulations that create video output to prevent mixing old frames with new ones.
Behavior
- Checks if the "video" directory exists
- Scans the directory for files with ".png" extension
- Attempts to delete each PNG file found
- Reports the number of files deleted
- Issues warnings for any files that cannot be deleted
Returns
Nothing
Example
# Clean up before creating new video frames
cleanup_video_folder()
# Run simulation that generates PNG frames
# ...
# Create video from frames
createVideo()
See Also
createVideo
: Create MP4 video from PNG framescreateAllVideos
: Create videos for all measurement types
FLORIDyn.createVideo
— FunctioncreateVideo(prefix::String; video_dir="video", output_dir="video", fps=2, delete_frames=false)
Convert PNG files in a directory starting with a given prefix into an MP4 video.
Arguments
prefix::String
: The prefix string that PNG files must start with (e.g., "velocityreduction", "windspeed")video_dir::String
: Directory containing the PNG files (default: "video")output_dir::String
: Directory where the output video will be saved (default: "video")fps::Int
: Frames per second for the output video (default: 2)delete_frames::Bool
: Whether to delete the PNG files after creating the video (default: false)
Returns
String
: Path to the created video file, or empty string if creation failed
Description
This function searches for PNG files in the specified directory that start with the given prefix, sorts them naturally (handling numeric sequences correctly), and combines them into an MP4 video using FFmpeg. The function requires FFmpeg to be installed on the system.
Examples
# Create video from velocity reduction frames
video_path = createVideo("velocity_reduction"; fps=4)
# Create video from wind speed frames and delete source frames
video_path = createVideo("wind_speed"; fps=6, delete_frames=true)
# Create video from custom directory
video_path = createVideo("added_turbulence"; video_dir="custom_plots", output_dir="videos")
Requirements
- FFmpeg must be installed and available in the system PATH
- PNG files should follow a consistent naming pattern with the prefix
- Recommended naming: "prefixt0000s.png", "prefixt0012s.png", etc.
Notes
- Files are sorted naturally to handle numeric sequences correctly (e.g., t0001s, t0010s, t0100s)
- The output video filename will be "prefix_animation.mp4"
- If no matching files are found, the function returns an empty string
- FFmpeg parameters are optimized for good quality and reasonable file size
FLORIDyn.createAllVideos
— FunctioncreateAllVideos(; video_dir="video", output_dir="video", fps=2, delete_frames=false)
Create videos for all common measurement types found in the video directory.
Arguments
video_dir::String
: Directory containing the PNG files (default: "video")output_dir::String
: Directory where output videos will be saved (default: "video")fps::Int
: Frames per second for output videos (default: 2)delete_frames::Bool
: Whether to delete PNG files after creating videos (default: false)
Returns
Vector{String}
: Paths to created video files
Description
This convenience function automatically detects common measurement type prefixes in the video directory and creates videos for each type found. It looks for the following prefixes:
- "velocity_reduction"
- "added_turbulence"
- "wind_speed"
Example
# Create videos for all measurement types found
video_paths = createAllVideos(fps=4, delete_frames=true)
println("Created videos: ", video_paths)
FLORIDyn.natural_sort_key
— Functionnatural_sort_key(filename::String)
Generate a sort key for natural sorting of filenames containing numbers.
Arguments
filename::String
: The filename to generate a sort key for
Returns
Vector
: Sort key that handles numeric sequences naturally
Description
This function creates a sort key that handles numeric sequences in filenames correctly. For example, it will sort ["file1.png", "file10.png", "file2.png"] as ["file1.png", "file2.png", "file10.png"] rather than alphabetically.
Examples
files = ["velocity_reduction_t0001s.png", "velocity_reduction_t0010s.png", "velocity_reduction_t0002s.png"]
sorted_files = sort(files, by=natural_sort_key)
# Result: ["velocity_reduction_t0001s.png", "velocity_reduction_t0002s.png", "velocity_reduction_t0010s.png"]