Low-Level Functions

Calculating the wind directions

FLORIDyn.getWindDirTFunction
getWindDirT(::Direction_Constant, wind_dir, iT, _)

Return wind direction in SOWFA-degrees for the requested turbine(s).

Arguments

  • wind_dir: The wind direction (scalar).
  • iT: Index or indices of the turbines (can be an integer or vector).
  • _: Placeholder for unused argument.

Returns

  • phi: Array of wind direction values, same size as iT.
source
getWindDirT(::Direction_Constant_wErrorCov, wind_dir::WindDirType, iT, t)

Return wind direction in SOWFA-deg for the requested turbine(s).

Arguments

  • wind_dir::WindDirType: WindDirType
  • iT: Vector of turbine indices (can be any indexable collection)
  • t: Time step

Returns

  • phi: Vector of wind directions for the selected turbines, including random perturbation
source
getWindDirT(::Direction_Interpolation, wind_dir::AbstractMatrix, iT, t)

Direction_Interpolation

Returns the wind direction at the respective turbine(s). Uniform interpolation version - all turbines experience the same changes.

Arguments:

  • wind_dir::Matrix: columns are time and phi (wind direction)
  • iT: single value or vector with turbine index/indices
  • t: time of request

Returns:

  • phi: Vector of wind directions for each turbine in iT [°]
source
getWindDirT(::Direction_Interpolation_wErrorCov, wind_dir::WindDirMatrix, iT, t)

Returns the wind direction at the respective turbine(s). Uniform interpolation version - all turbines experience the same changes.

Arguments:

  • wind_dir::WindDirMatrix: WindDirMatrix
  • iT: single value or vector with turbine index/indices
  • t: time of request

Returns:

  • phi: Vector of wind directions for each turbine in iT [°]
source
getWindDirT(::Direction_InterpTurbine, wind_dir, iT, t)

Return wind direction in SOWFA-degrees for the requested turbine(s).

Arguments

  • wind_dir::Matrix: Each row is [time, phi_T0, phi_T1, ...].
  • `iT: Index or indices of turbines.
  • t: Time of request. [s]

Returns

  • phi::Vector{Float64}: Wind direction(s) for the selected turbine(s) at time t. [°]
source
getWindDirT(::Direction_InterpTurbine_wErrorCov, wind_dir::WindDirMatrix, iT, t)

Return wind direction in SOWFA-deg for the requested turbine(s).

Arguments

  • wind_dir::WindDirMatrix: See: WindDirMatrix
  • iT: Index or indices of the turbines (can be integer or vector)
  • t: Time of request (Float64) [s]

Returns

  • phi: Wind direction(s) for requested turbine(s), perturbed with noise. [°]
source
getWindDirT(::Direction_RW_with_Mean, wind_dir_now, wind_dir::WindDirTriple)

Returns the wind direction at the respective turbine(s).

Arguments

  • wind_dir_now: Current value (vector)
  • wind_dir::WindDirTriple: WindDirTriple

Returns

  • phi: Updated wind direction(s) (vector) [°]
source
getWindDirT(::Direction_RW_with_Mean, wind_dir::WindDirTriple, iT, t)

Random walk with mean reversion model for wind direction.

Arguments

  • ::Direction_RW_with_Mean: Direction mode indicator
  • wind_dir::WindDirTriple: Wind direction data containing Init, CholSig, and MeanPull
  • iT: Turbine index or indices
  • t: Time value (unused in this implementation) [s]

Returns

  • phi: Wind direction(s) for the requested turbine(s) [°]
source
FLORIDyn.getWindDirT_EnKFFunction
getWindDirT_EnKF(::Direction_EnKF_InterpTurbine, wind_dir::AbstractMatrix, iT, t)

DirectionEnKFInterpTurbine

Return wind direction in SOWFA-deg for the requested turbine(s).

Arguments

  • wind_dir::Matrix: Matrix where each row is [time, phi_T0, phi_T1, ... phi_Tn]
  • iT: Index or indices of the turbines (can be integer or vector)
  • t: Time of request (scalar)

Returns

  • phi: Wind direction(s) at time t for turbine(s) iT [°]
source
FLORIDyn.getDataDirFunction
getDataDir(set::Settings, wind::Wind, wf::WindFarm, t)

Retrieve wind direction data for all turbines at the current simulation time.

Arguments

  • set::Settings: Simulation settings containing the direction mode configuration
  • wind::Wind: Wind field data structure containing direction information and input type
  • wf::WindFarm: Wind farm object containing turbine states and configuration
  • t: Current simulation time for temporal interpolation

Returns

  • phi: Wind direction values (typically in radians) for all turbines at the specified time

Description

This function reads wind direction data and returns the current wind direction angle (phi) for all turbines in the wind farm. The function handles two different input modes:

  1. Random Walk with Mean mode (wind.input_dir == "RW_with_Mean"): Uses the current wind farm state from wf.States_WF[wf.StartI, 2] along with the wind direction data to compute direction.

  2. Standard temporal interpolation mode: Uses the wind direction data directly with temporal interpolation for all turbines at the specified simulation time.

The function dispatches to getWindDirT with appropriate parameters based on the input mode, ensuring consistent wind direction estimation across different modeling approaches.

Examples

# Get wind direction for all turbines at current simulation time
phi = getDataDir(settings, wind_data, wind_farm, 100.0)

# The returned phi contains direction values for all turbines
direction_turbine_1 = phi[1]

Notes

  • The function automatically handles different wind input modes through conditional logic
  • For random walk mode, uses existing wind farm state as reference
  • For standard mode, performs temporal interpolation across all turbines

See also

  • getWindDirT: Underlying function for wind direction temporal interpolation
  • Settings: Configuration structure containing direction mode settings
  • Wind: Wind field data structure containing direction information and input type
  • WindFarm: Wind farm configuration structure
source
FLORIDyn.correctDir!Function
correctDir!(::Direction_All, set::Settings, wf::WindFarm, wind::Wind, t)

Apply direction correction to all turbines in the wind farm using the Direction_All strategy.

Arguments

  • ::Direction_All: The direction correction strategy type that applies corrections to all turbines
  • set::Settings: Simulation settings containing the direction mode configuration
  • wf::WindFarm: Wind farm object containing turbine states and configuration (modified in-place)
  • wind::Wind: Wind field data structure containing direction information and input type
  • t: Current simulation time for temporal interpolation

Returns

  • nothing: This function modifies the wind farm state in-place and returns nothing

Description

This function applies a direction correction using the Direction_All strategy, which updates the wind direction for all turbines in the wind farm. The function performs the following operations:

  1. Data Retrieval: Calls getDataDir to obtain current wind direction data for all turbines
  2. State Update: Updates the wind direction in the wind farm state (wf.States_WF[:, 2])
  3. Observation Point Orientation: If the state matrix has 4 columns, also updates the observation point orientation (wf.States_WF[wf.StartI, 4]) to match the wind direction

The correction is applied uniformly to all turbines using the first direction value from the retrieved direction data.

Examples

# Apply direction correction to all turbines
correctDir!(Direction_All(), settings, wind_farm, wind_data, 100.0)

# The wind farm state is modified in-place
current_direction = wind_farm.States_WF[1, 2]  # Updated direction for first turbine

Notes

  • This function modifies the wind farm state in-place (indicated by the ! suffix)
  • All turbines receive the same direction correction value (phi[1])
  • The observation point orientation is only updated if the state matrix has 4 columns
  • Direction values are typically in radians following standard wind engineering conventions

See also

  • getDataDir: Function for retrieving wind direction data
  • Direction_All: Direction correction strategy type
  • Settings: Simulation settings structure
  • WindFarm: Wind farm configuration structure
  • Wind: Wind field data structure
source

Calculating the wind velocity

FLORIDyn.getDataVelFunction
    getDataVel(set::Settings, wind::Wind, wf::WindFarm, t, tmp_m, floris::Floris)

Return the wind speed vector u for all turbines at simulation time t according to the configured input / model mode. Also returns (potentially updated) wind.

Arguments

  • set::Settings: simulation settings (uses set.vel_mode)
  • wind::Wind: wind field state (wind.input_vel selects special branches)
  • wf::WindFarm: wind farm (uses wf.nT, wf.States_WF)
  • t: current simulation time
  • tmp_m: temporary matrix (only used for IandI wake reduction)
  • floris::Floris: FLORIS parameters (yaw exponent etc., only IandI branch)

Supported (unit tested in test_getDataVel_branches.jl)

  • Standard interpolation / constant variants via set.vel_mode: Velocity_Constant, Velocity_Interpolation, Velocity_Constant_wErrorCov, Velocity_Interpolation_wErrorCov, Velocity_InterpTurbine, Velocity_InterpTurbine_wErrorCov, Velocity_ZOH_wErrorCov.
  • EnKF turbine interpolation branch: wind.input_vel == "EnKF_InterpTurbine" calling getWindSpeedT_EnKF(Velocity_EnKF_InterpTurbine(), ...) with clamping of out-of-range times.

Not yet fully integrated (guarded / broken tests)

  • "I_and_I": an internal estimator state struct (WSEStruct in windfield_velocity.jl) already exists and the low-level update routine WindSpeedEstimatorIandI_FLORIDyn runs, but a public, documented construction path (export, convenience constructor, validation of required fields, tests) is missing. The branch is therefore kept experimental and the test remains @test_broken until we provide a stable API (e.g. build_IandI_estimator(wf, data; kwargs...)).
  • "RW_with_Mean": random-walk-with-mean model commented out; current call raises MethodError.

Planned cleanups / TODO

  1. Provide concrete exported estimator type & finalize I_and_I logic.
  2. Re-introduce Random Walk with Mean model (Velocity_RW_with_Mean).

Behavior summary

  • Default: u = getWindSpeedT(set.vel_mode, wind.vel, 1:nT, t).
  • EnKF: per-turbine linear interpolation table with time clamping.
  • IandI: (future) estimator integration plus optional wake reduction using tmp_m[:,1].
  • RWwithMean: (future) stochastic update around mean with mean-pull term.

All returned velocities are in m/s. Only IandI may mutate wind.vel estimator state.

Example

u, wind = getDataVel(set, wind, wf, 100.0, tmp_m, floris)
source

Calculating the wind shear

FLORIDyn.getWindShearTFunction
getWindShearT(::Shear_Interpolation, wind_shear::AbstractMatrix, z)

Compute the wind shear at a given height z using the specified wind_shear model.

Arguments

  • ::Shear_Interpolation: (Type only) Use interpolation to determine the wind shear.
  • wind_shear: A matrix describing the wind shear profile.
  • z: The height (in meters) at which to evaluate the wind shear.

Returns

  • The wind shear value at height z.

REMARKS

Expects a .csv file called "WindShearProfile.csv" with a normalized wind speed profile for different heights:

z, (u_z/u0)
z, (u_z/u0)
z, (u_z/u0)

There is a linear interpolation between every pair. In case z is out of bounds the function will use the closest available setpoint.

source
getWindShearT(::Shear_PowerLaw, wind_shear::WindShear, z_norm)

Return the shear factor u_eff = shear * u_referenceHeight using the power law.

Arguments

  • Shear_PowerLaw: (type only, unused) Specifies that this method applies to the power law model
  • wind_shear: A struct of type (WindShear)(@ref)
    • z0: Reference height (not used in this function)
    • alpha: WindShear coefficient
  • z_norm: Height(s) (can be scalar or array)

Returns

  • shear: The shear factor at the given height(s)
source

Calculating the wind turbulence

FLORIDyn.getDataTIFunction
getDataTI(set::Settings, wind::Wind, wf::WindFarm, t) -> Vector

Retrieve turbulence intensity data for all turbines at the current simulation time.

This function obtains turbulence intensity values for all turbines in the wind farm using the configured turbulence model and wind field data. It serves as a wrapper around the underlying turbulence intensity retrieval system.

Arguments

  • set::Settings: Settings object containing simulation configuration
    • set.turb_mode: Turbulence model configuration specifying the retrieval method
  • wind::Wind: Wind configuration object containing turbulence intensity data
    • wind.ti: Turbulence intensity data, parameters, or model configuration
  • wf::WindFarm: Wind farm object containing turbine information
    • wf.nT: Number of turbines in the wind farm
  • t: Current simulation time for time-dependent turbulence intensity models

Returns

  • Vector: Turbulence intensity values for all turbines (dimensionless, typically 0.05-0.25)

Behavior

The function creates a vector of all turbine indices [1, 2, ..., nT] and retrieves the corresponding turbulence intensity values using the specified turbulence model. The actual retrieval method depends on the set.turb_mode configuration and can include:

  • Constant turbulence intensity
  • Time-interpolated values from data files
  • Turbine-specific interpolation
  • Random walk models with covariance

Example

# Get turbulence intensity for all turbines at t=100s
TI_values = getDataTI(settings, wind_config, wind_farm, 100.0)
println("TI for turbine 1: ", TI_values[1])

See Also

  • getWindTiT: Underlying function for turbulence intensity retrieval
  • correctTI!: Function that uses this data to update wind farm states
source
FLORIDyn.getWindTiTFunction
getWindTiT(::TI_Constant, wind_ti, iT, _)

Return turbulence intensity for the requested turbine(s).

Arguments

  • ::TI_Constant: Type parameter to indicate constant wind turbulence
  • wind_ti: Constant value (turbulence intensity)
  • iT: Index or indices of the turbines
  • _: will be ignored

Returns

  • Ti: Array of turbulence intensity values for each turbine index
source
getWindTiT(::TI_Interpolation, wind_ti::AbstractMatrix, iT, t)

Interpolates the wind turbulence intensity (TI) at a given time t using the specified TI_Interpolation method.

Arguments

  • ::TI_Interpolation: Use linear interpolation to calculate the turbulence intensity.
  • wind_ti::Matrix: Matrix containing wind turbulence intensity values over time.
  • iT: Index/indices of the turbines (can be Int or array).
  • t: The specific time at which to interpolate the turbulence intensity.

Returns

  • The interpolated turbulence for the requested turbine(s) at time t.

Notes

  • The function assumes that wind_ti contains the necessary data for interpolation as (time, TI) pairs (n×2 matrix)
  • Uniform interpolation version - all turbines experience the same changes.
source
getWindTiT(::TI_InterpTurbine, wind_ti::AbstractMatrix, iT, t)

Retrieve the wind turbulence intensity (TI) for a specific turbine at a given time.

Arguments

  • ::TI_InterpTurbine: The turbulence intensity interpolation object for the turbine.
  • wind_ti::AbstractMatrix: Matrix containing wind turbulence intensity values.
  • iT: Index of the turbine for which the TI is requested.
  • t: Time at which the TI value is needed.

Returns

  • The interpolated wind turbulence intensity value for the specified turbine at time t.
source
FLORIDyn.correctTI!Function
correctTI!(::TI_None, set::Settings, wf::WindFarm, wind::Wind, t) -> Nothing

Update turbulence intensity values in the wind farm state matrix without correction.

This function implements the "no correction" strategy for turbulence intensity, where the wind farm turbulence intensity values are updated with fresh data from the wind field model without applying any correction algorithms. It serves as the baseline approach for turbulence intensity handling in FLORIDyn simulations.

Arguments

  • ::TI_None: Dispatch type indicating no turbulence intensity correction algorithm
  • set::Settings: Settings object containing simulation configuration and turbulence model parameters
    • set.turb_mode: Turbulence model configuration specifying the retrieval method
  • wf::WindFarm: Wind farm object containing the state matrices to be updated
    • wf.States_WF: Wind field states matrix where column 3 contains turbulence intensity values
    • wf.StartI: Starting indices for each turbine's observation points
    • wf.nT: Number of turbines
  • wind::Wind: Wind configuration object containing turbulence intensity data
    • wind.ti: Turbulence intensity data or model parameters
  • t: Current simulation time for time-dependent turbulence intensity retrieval

Returns

  • Nothing: The function modifies the wind farm state in-place

Behavior

  1. Retrieves current turbulence intensity values for all turbines using getDataTI
  2. Updates the wind farm state matrix wf.States_WF at rows wf.StartI and column 3
  3. Transposes the turbulence intensity vector to match the matrix structure
  4. Provides error handling for matrix update operations

Example

# Update turbulence intensity without correction at t=50s
correctTI!(TI_None(), settings, wf, wind, 50.0)

# The wind farm state matrix is now updated with new TI values
current_ti = wf.States_WF[wf.StartI, 3]

Notes

  • The function modifies the wind farm object in-place (indicated by the ! suffix)
  • This "no correction" approach provides baseline turbulence intensity without applying wake-induced corrections or measurement-based adjustments
  • Error handling ensures graceful failure if matrix dimensions are incompatible

See Also

  • getDataTI: Function used to retrieve turbulence intensity data
  • TI_None: Dispatch type for no correction strategy
source

Controller functions

FLORIDyn.getYawFunction
getYaw(::Yaw_SOWFA, con_yaw_data::AbstractMatrix, iT, t) -> Float64 or Vector{Float64}

Return the yaw angle at time t for the specified turbine(s) using linear interpolation.

Arguments

  • ::Yaw_SOWFA: Controller type dispatch parameter for SOWFA-style yaw control
  • con_yaw_data::Matrix{Float64}: Control data matrix where:
    • First column contains time values (in seconds)
    • Subsequent columns contain yaw angles for each turbine (in degrees)
  • iT: Turbine index or indices to query:
    • Integer: Single turbine index (1-based)
    • AbstractVector{<:Integer}: Vector of turbine indices for multiple turbines
  • t::Real: Requested time (in seconds)

Returns

  • Float64: Single yaw angle (in degrees) if iT is an integer
  • Vector{Float64}: Vector of yaw angles (in degrees) if iT is a vector

Behavior

  • Interpolation: Uses linear interpolation between time points with flat extrapolation
  • Out-of-bounds handling: If t is outside the time range, the function:
    • Issues a warning message
    • Clamps t to the nearest boundary (first or last time point)
  • Single time point: If only one time point exists, returns the corresponding yaw value directly
  • Error handling: Throws an error if iT is not an integer or vector of integers

Data Format

The con_yaw_data matrix should have the structure:

[time₁  yaw₁₁  yaw₁₂  ...  yaw₁ₙ]
[time₂  yaw₂₁  yaw₂₂  ...  yaw₂ₙ]
[  ⋮      ⋮      ⋮    ⋱     ⋮  ]
[timeₘ  yawₘ₁  yawₘ₂  ...  yawₘₙ]

where m is the number of time steps and n is the number of turbines.

Examples

# Example control data: 3 time points, 2 turbines
con_yaw_data = [0.0  10.0  5.0;   # t=0s: T1=10°, T2=5°
                1.0  15.0  10.0;  # t=1s: T1=15°, T2=10°
                2.0  20.0  15.0]  # t=2s: T1=20°, T2=15°

# Get yaw for turbine 1 at t=0.5s (interpolated)
yaw1 = getYaw(Yaw_SOWFA(), con_yaw_data, 1, 0.5)  # Returns 12.5°

# Get yaw for multiple turbines at t=1.5s
yaws = getYaw(Yaw_SOWFA(), con_yaw_data, [1, 2], 1.5)  # Returns [17.5°, 12.5°]

# Out-of-bounds time (will issue warning)
yaw_oob = getYaw(Yaw_SOWFA(), con_yaw_data, 1, 5.0)  # Returns 20.0° with warning

See Also

  • Yaw_SOWFA: Controller type for SOWFA-style yaw control
  • Interpolations.linear_interpolation: Underlying interpolation method used
source
getYaw(::Yaw_Constant, con_yaw_data::AbstractMatrix, iT, t)

Return a single constant yaw angle (degrees) for one or multiple turbines.

Assumptions / Simplified Layout:

  • con_yaw_data is a matrix with at least one row and one column.
  • Only the value con_yaw_data[1,1] is used; any extra rows/columns are ignored.

Arguments:

  • ::Yaw_Constant: dispatch marker
  • con_yaw_data: matrix whose first element holds the constant yaw angle
  • iT: turbine index (Integer) or vector of indices
  • t: ignored (kept for a uniform interface)

Returns:

  • Float64 when iT is an Integer
  • Vector{Float64} when iT is a vector
source

Helper functions

FLORIDyn.toMSRFunction
toMSR(s::String)

Converts the input string s to a MSR (Measurement System Representation) enumeration.

Supports the following string formats:

  • Enum names: "VelReduction", "AddedTurbulence", "EffWind"
  • Flow field names: "flow_field_vel_reduction", "flow_field_added_turbulence", "flow_field_eff_wind_speed"
  • Measurement names: "msr_vel_reduction", "msr_added_turbulence", "msr_eff_wind_speed"

See also MSR.

source
FLORIDyn.now_microsecondsFunction
now_microseconds()::String

Returns current timestamp as a string with microsecond resolution. Format: "YYYY-mm-ddTHH-MM-SS.uuuuuu"

Examples

julia> now_microseconds()
"2025-08-08T16-58-55.494911"
source
FLORIDyn.now_nanosecondsFunction
now_nanoseconds()::String

Returns current timestamp as a string with nanosecond resolution. Format: "YYYY-mm-ddTHH-MM-SS.nnnnnnnnn"

Examples

julia> now_nanoseconds()
"2025-08-08T16-58-55.494911123"
source
FLORIDyn.precise_nowFunction
precise_now()::String

Alias for now_microseconds() - returns current timestamp with microsecond resolution.

source
FLORIDyn.unique_nameFunction
unique_name()::String

Creates a unique directory name for storing simulation results.

This function generates a unique directory name by combining the prefix "floridynrun" with a high-resolution timestamp. The timestamp includes microsecond precision to ensure uniqueness even when multiple simulations are started in quick succession.

Returns

  • String: A unique directory name in the format "floridyn_run_YYYY-mm-ddTHH-MM-SS.uuuuuu"

Examples

julia> unique_name()
"floridyn_run_2025-08-08T16-58-55.494911"

julia> unique_name()
"floridyn_run_2025-08-08T16-58-55.495123"

Notes

  • The generated name is suitable for use as a directory name on all operating systems
  • Uses hyphens instead of colons in the time portion for filesystem compatibility
  • Microsecond precision ensures uniqueness for rapid successive calls
  • Used to create separate output directories for each run

See Also

  • now_microseconds: The underlying timestamp function used
  • Vis: Visualization settings that may use unique names for output directories
source
FLORIDyn.isdelftblueFunction
isdelftblue() -> Bool

Check if the current environment is the Delft Blue supercomputer.

This function determines whether the code is running on the Delft Blue supercomputer by checking for the existence of the ~/scratch directory, which is a characteristic feature of the Delft Blue file system.

Returns

  • Bool: true if running on Delft Blue (i.e., ~/scratch directory exists), false otherwise.

Usage

This function is used throughout FLORIDyn.jl to automatically adapt file paths to the computing environment. On Delft Blue, output files are stored in the scratch directory, which can be accessed remotely and offers lots of space.

Example

if isdelftblue()
    output_path = joinpath(homedir(), "scratch", "out")
else
    output_path = joinpath(pwd(), "out")
end
source
FLORIDyn.find_floridyn_runsFunction
find_floridyn_runs(directory::String = pwd())

Find all directories starting with "floridyn_run" in the specified directory.

Arguments

  • directory::String: Directory to search in (default: current working directory)

Returns

  • Vector{String}: List of full paths to floridyn_run directories, sorted by name

Examples

# Find floridyn_run directories in current directory
dirs = find_floridyn_runs()

# Find in specific directory
dirs = find_floridyn_runs("out")

# Count how many exist
count = length(find_floridyn_runs())

See Also

source
FLORIDyn.delete_resultsFunction
delete_results(vis::Vis, n::Int=1; dry_run::Bool = false)

Delete the newest n directories starting with "floridyn_run" from both the visualization output directory and video directory.

This function provides comprehensive cleanup by removing the most recent floridyn_run directories from both output and video paths simultaneously. It's particularly useful for removing failed runs, test runs, or managing disk space by keeping only the most relevant simulation outputs across both directory types.

Arguments

  • vis::Vis: Visualization settings object containing both output and video path configurations
  • n::Int: Number of newest directories to delete from each directory (must be positive, default: 1)
  • dry_run::Bool: Preview mode - shows what would be deleted without actually deleting (default: false)

Returns

  • Vector{String}: Absolute paths of directories that were deleted from both locations combined. Returns empty vector if no matching directories found or if n ≤ 0.

Behavior

  • Directory Search: Searches both vis.output_path and vis.video_path for directories matching floridyn_run_* pattern
  • Dual Cleanup: Operates on both output and video directories in sequence
  • Sorting: Sorts directories by modification time (newest first) within each directory
  • Selection: Selects up to n newest directories for deletion from each location
  • Dry Run: When dry_run=true, logs what would be deleted but performs no actual deletion

Examples

# Create visualization settings
vis = Vis("data/vis_default.yaml")

# Delete the single newest floridyn_run directory from both output and video directories
deleted = delete_results(vis)
println("Deleted: ", length(deleted), " directories total")

# Delete the 3 newest floridyn_run directories from each location
deleted = delete_results(vis, 3)
println("Deleted directories: ", basename.(deleted))

Important Notes

  • Dual Operation: This function operates on BOTH output and video directories
  • Independent Processing: Each directory is processed separately - n directories from output AND n directories from video
  • Non-existent Directories: Silently skips directories that don't exist rather than creating them

See Also

  • Vis: Visualization settings struct with output and video path configuration
  • unique_name(): Creates timestamped floridyn_run directories
  • find_floridyn_runs(): Lists existing floridyn_run directories in a given path
source