Reference

FLOWFarm.DiscretizedWindResourceType
DiscritizedWindResource(wind_directions, wind_speeds, wind_probabilities, measurement_heights, air_density, ti_model, wind_shear_model)

Struct defining a wind resource

Arguments

  • wind_directions::Array{Float,1}(Nstates): an array of wind directions corresponding to each wind farm state in radians
  • wind_speeds::Array{Float,1}(Nstates): an array of wind speeds corresponding to each wind farm state in meters/second
  • wind_probabilities::Array{Float,1}(Nstates): an array of probabilities corresponding to each wind farm state with values in [0,1]
  • measurement_heights::Array{Float,1}(Nstates): an array of measurement heights corresponding to each wind farm state
  • air_density::Float: the air density
  • ambient_ti::Array{Float,1}: an array of the ambient turbulence intensity for each wind direction
  • wind_shear_model::Array{AbstractWindShearModel}(1): contains a struct defining the desired turbulence intensity model
source
FLOWFarm.GaussOriginalType
GaussOriginal(k_star)

Container for parameters related to the origina Gaussian deficit model presented by Bastankhah and Porte-Agel 2014

Arguments

  • k_star::Float: parameter controlling the wake spreading rate and deficit decay. Default value is 0.075
source
FLOWFarm.GaussSimpleType
GaussSimple(k, wec_factor)

Container for parameters related to the Gaussian deficit model with yaw presented by Bastankhah and Porte-Agel 2016

Arguments

  • k::Float: parameter controlling the spread of the wake
  • wec_factor::Array{Float}: parameter artificial wake spreading for wake expansion continuation (WEC) optimization
source
FLOWFarm.GaussYawType
GaussYaw(turbulence_intensity, horizontal_spread_rate, vertical_spread_rate, alpha_star, beta_star)

Container for parameters related to the Gaussian deficit model with yaw presented by Bastankhah and Porte-Agel 2016

Arguments

  • horizontal_spread_rate::Float: parameter controlling the horizontal spread of the deficit model. Default value is 0.022.
  • vertical_spread_rate::Float: parameter controlling the vertical spread of the deficit model. Default value is 0.022.
  • alpha_star::Float: parameter controlling the impact of turbulence intensity on the length of the near wake. Default value is 2.32.
  • beta_star::Float: parameter controlling the impact of the thrust coefficient on the length of the near wake. Default value is 0.154.
source
FLOWFarm.GaussYawDeflectionType
GaussYawDeflection(horizontal_spread_rate, vertical_spread_rate, alpha_star, beta_star)

Container for parameters related to the Gaussian deflection model presented by Bastankhah and Porte-Agel 2016

Arguments

  • horizontal_spread_rate::Float: parameter controlling the horizontal spread of the deficit model. Default value is 0.022.
  • vertical_spread_rate::Float: parameter controlling the vertical spread of the deficit model. Default value is 0.022.
  • alpha_star::Float: parameter controlling the impact of turbulence intensity on the length of the near wake. Default value is 2.32.
  • beta_star::Float: parameter controlling the impact of the thrust coefficient on the length of the near wake. Default value is 0.154.
source
FLOWFarm.GaussYawVariableSpreadType
GaussYawVariableSpread(turbulence_intensity, horizontal_spread_rate, vertical_spread_rate, alpha_star, beta_star)

Container for parameters related to the Gaussian deficit model with yaw presented by Bastankhah and Porte-Agel 2016 and the farm model presented by Niayifar and Porte-Agel in 2016.

Arguments

  • alpha_star::Float: parameter controlling the impact of turbulence intensity on the length of the near wake. Default value is 2.32.
  • beta_star::Float: parameter controlling the impact of the thrust coefficient on the length of the near wake. Default value is 0.154.
source
FLOWFarm.GaussYawVariableSpreadDeflectionType
GaussYawDeflectionVariableSpread(alpha_star, beta_star, k1, k2, wec_factor)

Container for parameters related to the Gaussian deflection model with yaw presented by Bastankhah and Porte-Agel 2016

Arguments

  • alpha_star::Float: parameter controlling the impact of turbulence intensity on the length of the near wake. Default value is 2.32.
  • beta_star::Float: parameter controlling the impact of the thrust coefficient on the length of the near wake. Default value is 0.154.
  • k1::Float: first parameter tuning wake spread as based on turbulence intensity
  • k2::Float: second parameter tuning wake spread as based on turbulence intensity
source
FLOWFarm.JensenCosineType
JensenCosine(alpha)

Container for parameters related to the Jensen Cosine deficit model

Arguments

  • alpha::Float: parameter controlling the wake deficit decay rate. Default value is 0.1
  • beta::Float: parameter controlling the width of the cosine function. Default value is 20.0 deg., given in radians.
source
FLOWFarm.JensenTopHatType
JensenTopHat(alpha)

Container for parameters related to the Jensen Top Hat deficit model

Arguments

  • alpha::Float: parameter controlling the wake spreading rate and deficit decay. Default value is 0.1
source
FLOWFarm.JiminezYawDeflectionType
JiminezYawDeflection(horizontal_spread_rate)

Container for parameters related to the Jiminez deflection model

Arguments

  • horizontal_spread_rate::Float: parameter controlling the wake spreading rate and deficit decay. Default value is 0.1
source
FLOWFarm.LevelizedType
Levelized(TCC, BOS, FC, FCR, OpEx)

Container for parameters related to the Levelized Cost of Energy model (NREL 2016 Cost of Wind Energy)

Arguments

Arguments

  • TCC::Float: Turbine Capital Cost not including the tower module
  • BOS::Float: Balance of System (Costs outside of turbine i.e. operation and maintenance)
  • FC::Float: Financial Costs including construction and contingency
  • FCR::Float: Fixed Charge Rate
  • OpEx::Float: Operational Expenditures
source
FLOWFarm.LocalTIModelMaxTIType
LocalTIModelMaxTI(astar, bstar, k1, k2)

Calculate local turbulence intensity using the model presented in Niayifar and Porte Agel (2015, 2016)

Arguments

  • astar::Float: wake spreading parameter from Bastankhah and Porte-Agel Gaussian wake model
  • bstar::Float: wake spreading parameter from Bastankhah and Porte-Agel Gaussian wake model
  • k1::Float: slope of k vs TI curve
  • k2::Float: vertical offset of k vs TI curve
source
FLOWFarm.MultiZoneType
MultiZone(me, ke, MU, aU, bU)

Container for parameters related to the MultiZone deficit model

Arguments

  • me::Float: parameter controlling general wake expansion. Default value is 0.065
  • ke::Array{Float}(3): parameters controlling the wake expansion of each zone respectively. Default values are [-0.5 0.22 1.0].
  • MU::Array{Float}(3): parameters controlling the wake deficit decay of each zone respectively. Default values are [0.5 1.0 5.5].
  • aU::Float: parameter impacting the wake deficit decay for a constant wake deflection. Default value is 5.0.
  • bU::Float: parameter changing the wake deficit decay under yawed conditions. Default value is 1.66.
source
FLOWFarm.MultizoneDeflectionType
MultizoneDeflection(horizontal_spread_rate)

Container for parameters related to the Jiminez deflection model

Arguments

  • horizontal_spread_rate::Float: parameter controlling the wake spreading rate and deficit decay. Default value is 0.1
  • ad::Float:Helps define the horizontal deflection of the wake at 0 deg yaw
  • bd::Float:Helps define the horizontal deflection of the wake due to downwind distance at 0 deg yaw
source
FLOWFarm.PowerLawWindShearType
PowerLawWindShear(shear_exponent, ground_height)

Provides shear exponent and ground height to define wind shear curve. Ground height may be tuned because the power law does not always hold near the ground.

Arguments

  • shear_exponent::Float: defines trajectory of wind shear
  • ground_height::Float: height of the ground (typically zero)
  • shear_order::Bool: when shear should be calculated. Can be "first", "last", or "none"
source
FLOWFarm.PowerModelConstantCpType
PowerModelConstantCp(cp)

Models will assume a constant cp value as provided

Arguments

  • cp::Float: constant power coefficient value
  • 'pp::TI': exponent for adjusting power for wind turbine yaw
source
FLOWFarm.PowerModelCpPointsType
PowerModelCpPoints(vel_points, cp_points)

Models will use adjust cp based on cp curve using linear interpolation of provided points

Arguments

  • vel_points::Array{N,Float}: wind speed values in m/s
  • cp_points::Array{N,Float}: power coefficient values corresponding to the provided speeds
  • 'pp::TF': exponent for adjusting power for wind turbine yaw
source
FLOWFarm.PowerModelPowerCurveCubicType
PowerModelPowerCurveCubic()

Power will be calculated based on turbine specifications assuming a cubic power curve. Note that this method is inherently incorrect and should only be used for theoretical purposes or after careful validation.

Arguments

  • 'pp::TF': exponent for adjusting power for wind turbine yaw
source
FLOWFarm.PowerModelPowerPointsType
PowerModelPowerPoints(vel_points, cp_points)

Models will use adjust wind turbine power based on power curve using linear interpolation of provided points

Arguments

  • vel_points::Array{N,Float}: wind speed values in m/s
  • power_points::Array{N,Float}: power values corresponding to the provided speeds
  • 'pp::TF': exponent for adjusting power for wind turbine yaw
source
FLOWFarm.ThrustModelConstantCtType
ThrustModelConstantCt(ct::Float)

Stores a constant ct value for wake calculations

Arguments

  • ct::Float: a constant ct value for computation
source
FLOWFarm.ThrustModelCtPointsType
ThrustModelCtPoints(vel_points, ct_points)

Stores the thrust coefficient curve in terms of corresponding velocity and thrust coefficient points. ct and velocity points should be in the same order and ordered from lowest wind speed to highest wind speed.

Arguments

  • inflow_velocity::Float: inflow velocity of the wind turbine
  • thrust_model::ThrustModelCtPoints: Struct containing ct and velocity points for ct curve
source
FLOWFarm.WindFarmType

WindFarm(windfarm, windresource, windfarmstates)

Struct defining a wind farm

Arguments

  • turbine_x::Array{Float}(Nturbines): contains windturbine x coordinates in the global reference frame
  • turbine_y::Array{Float}(Nturbines): contains windturbine y coordinates in the global reference frame
  • turbine_z::Array{Float}(Nturbines): contains windturbine base/z coordinates in the global reference frame
  • turbine_definition_ids::Array{Int}(Nturbines): contains integers for each wind turbine specifying its definition
  • turbine_definitions::Array{AbstractTurbineDefinition}(Ntypes): contains structs defining each wind turbine definition (design) used in the farm
source
FLOWFarm.WindFarmModelSetType
WindFarmModelSet(wakedeficitmodel, wake_deflection_model, wake_combination_model, local_ti_model)

Container for objects defining models to use in wind farm calculations

Arguments

  • wake_defiict_model::AbstractWakeDeficitModel: contains a struct defining the desired wake deficit model
  • wake_deflection_model::AbstractWakeDeflectionModel: contains a struct defining the desired wake deflection model
  • wake_combination_model::AbstractWakeCombinationModel: contains a struct defining the desired wake combination model
  • local_ti_model::AbstractTurbulenceIntensityModel: contains a struct defining the desired turbulence intensity model
source
FLOWFarm.DiscreteCircumFunction
PointsOnCircum(center_x, center_y, r, n = 100)

Given a circle center, radius, and number of discrete points, returns an array of discrete points along the circle's circumference

Arguments

  • center_x::Float64 : cartesian x-coordinate for the center of the circle
  • center_y::Float64 : cartesian y-coordinate for the center of the circle
  • r::Float64 : distance from circle's center to the circumference points
  • n::Float64 : defaults to 100, is the number of discrete evenly-spaced points that will be returned along the circle's circumference
source
FLOWFarm.GaussianTIMethod
GaussianTI(loc,turbine_x, turbine_y, rotor_diameter, hub_height, turbine_ct, 
    sorted_turbine_index, ambient_ti; div_sigma=2.5, div_ti=1.2)

Calculate local turbulence intensity based on "On wake modeling, wind-farm gradients and AEP predictions at the Anholt wind farm" by Pena Diaz, Alfredo; Hansen, Kurt Schaldemose; Ott, Søren; van der Laan, Paul ??

Arguments

  • loc::Array{Float,3}: [x,y,z] location of point of interest in wind direction ref. frame
  • turbine_x::Array{Float,nTurbines}: turbine wind direction locations in the wind direction reference frame
  • turbine_y::Array{Float,nTurbines}: turbine cross wind locations in the wind direction reference frame
  • rotor_diameter::Array{Float,nTurbines}: rotor diameters of all turbines
  • hub_height::Array{Float,nTurbines}: hub heights of all turbines relative to the ground
  • turbine_ct::Array{Float,nTurbines}: thrust coefficient of each turbine for the given state
  • sorted_turbine_index::Array{Float,nTurbines}: turbine north-south locations in the global reference frame
  • ambient_ti::Float: ambient turbulence intensity
  • div_sigma::Float: ?
  • div_ti::Float: ?
source
FLOWFarm.VR_boundaryMethod
VR_boundary(bndry_x_clsd, bndry_y_clsd, start_dist, turb_spacing, num_turbs, bndry_seg_length)

Uses the Boundary portion of Boundary-Grid variable reduction method place turbines along a closed wind farm boundary and perturb their location with one (1) variable <startdist>. NOTE: Use of this function assumes prior use of VRbounary_startup(), which ensures the number of turbines placed on the boundary doesn't violate any minimum spacing rules eiter along the boundary or around corners.

Arguments

  • bndry_x::Array{Float,1} : 1-D array of x-coordinates for the vertices around a singlar closed boundary
  • bndry_y::Array{Float,1} : 1-D array of y-coordinates for the vertices around a singlar closed boundary
  • start_dist::Float64: the distance (positive or negative) along the boundary from the first boundary point where the turbines will begin to be placed
  • turb_spacing::Float64: the fixed distance along the boundary's edge between adjacent turbines
  • 'numturbs::Float64`: the number of turbines to be placed around the boundary. Note that this function assumes VRbounary_startup() has already been run so that the user won't attempt to place too many turbines.
source
FLOWFarm.VR_boundary_startupMethod
VR_bounary_startup(bndry_x_clsd, bndry_y_clsd, start_dist, turb_min_spacing, num_turbs)

Determines if the requested number of turbines can be placed along the closed boundary with spacing and corner constraints. If the requested <num_turbs> is too many, places as many turbines as possible along the boundary, and returns the number of turbines not placed. NOTE: A shortcoming is that the smallest-angled corner limits the spacing of all turbines. in the worst case, a very thin boundary area would prevent any more than one turbine being placed on the boundary, though more would be optimal. Future work would check to make sure this corner (and the length of its adjacent sides) don't actually require limiting the minimum distance between turbines.

Arguments

  • bndry_x::Array{Float,1} : 1-D array of x-coordinates for the vertices around a singlar closed boundary
  • bndry_y::Array{Float,1} : 1-D array of y-coordinates for the vertices around a singlar closed boundary
  • start_dist::Float64: the distance (positive or negative) along the boundary from the first boundary point where the turbines will begin to be placed
  • turb_min_spacing::Float64: the fixed distance along the boundary's edge between adjacent turbines
  • 'numturbs::Float64`: the number of turbines to be placed around the boundary. Note that this function assumes VRbounary_startup() has already been run so that the user won't attempt to place too many turbines.
  • 'bndryseglength::Array{Int}`: an array of the lengths between adjacent boundary verticies, corresponding to how they appear in bndry_x and _y
source
FLOWFarm._ct_to_axial_ind_funcMethod
_ct_to_axial_ind_func(ct)

Calculate axial induction from the thrust coefficient. See Gebraad et. al. 2017 "Maximization of the Annual Energy Production of Wind Power Plants by Optimization of Layout and Yaw-Based Wake Control"

Arguments

  • ct::Float: thrust coefficient
source
FLOWFarm._gauss_yaw_potential_coreMethod
_gauss_yaw_potential_core(dt, yaw, ct, as, ti, bs)

Helper function for wakedeficitmodel when using the GaussYaw model. Computes the length of the near wake potential core.

source
FLOWFarm._gauss_yaw_spreadMethod
_gauss_yaw_spread(dt, k, dx, x0, yaw)

Helper function for wakedeficitmodel when using the GaussYaw model. Computes the standard deviation of the wake.

source
FLOWFarm._gauss_yaw_spread_interpolatedMethod
_gauss_yaw_spread_interpolated(dt, k, dx, x0, yaw)

Helper function for wakedeficitmodel when using the GaussYaw model. Computes the standard deviation of the wake. with an interpolation on the near wake.

source
FLOWFarm._niayifar_added_ti_functionMethod
_niayifar_added_ti_function(x, d_dst, d_ust, h_ust, h_dst, ct_ust, kstar_ust, delta_y, 
    ti_amb, ti_ust, ti_dst, ti_area_ratio_in; s=700.0)

Main code for calculating the local turbulence intensity at a turbine using the method of Niayifar and Porte Agel (2015, 2016).

Arguments

  • x::Float: downstream distance from turbine to point of interest
  • d_dst::Float: downstream turbine rotor diameter
  • d_ust::Float: upstream turbine rotor diameter
  • h_ust::Float: upstream turbine hub height
  • h_dst::Float: downstream turbine hub height
  • ct_ust::Float: upstream turbine thrust coefficient
  • kstar_ust::Float: upstream turbine wake expansion rate
  • delta_y::Float: cross wind separation from turbine to point of interest
  • ti_amb::Float: ambient turbulence intensity
  • ti_ust::Float: upstream turbine local turbulence intensity
  • ti_dst::Float: downstream turbine local turbulence intensity
  • ti_area_ratio_in::Float: current value of TI-area ratio for use in calculatin local TI
  • s::Float: smooth max smootheness parameter
source
FLOWFarm._remove_out_of_bounds_pointsMethod
_remove_perimeter_points!(n; alpha=0.0)

Internal function. Removes points outside or outside and on the border of the rotor-swept area

Arguments

  • y::AbstractArray: horizontal point locations
  • z::AbstractArray: vertical point locations
  • use_perimeter_points::Bool: flag that determines whether or not to include points on the boundary of the rotor-swept area
source
FLOWFarm.add_turbine!Method
add_turbine!(ax; view="side", hubdiameter=0.1, hubheight=0.9, radius=0.5, chord=0.1, 
    nacellewidth=0.3, nacelleheight=0.1, towerbottomdiam=0.1, towertopdiam=0.05, 
    overhang=0.05, s=5)

Convenience function for adding wind turbines to plots.

Arguments

  • ax::PyCall.PyObject: pre-initialized axis from pyplot
  • view::Number: determines which turbine view to use "top" or "side" (default)
  • hubdiameter::Number: hub diameter in axis coordinate frame
  • hubheight::Number: hub height in axis coordinate frame
  • radius::Number: full rotor radius in axis coordinate frame
  • chord::Number: maximum chord in axis coordinate frame
  • nacellewidth::Number: nacelle width in axis coordinate frame
  • nacelleheight::Number: nacelle height in axis coordinate frame
  • towerbottomdiam::Number: tower bottom diameter in axis coordinate frame
  • towertopdiam::Number: tower top diameter in axis coordinate frame
  • overhang::Number: overhang (distance from blade attachment to tower bottom in x axis) in axis coordinate frame
  • s::Number: scales overhang and tower location in x direction to work with condensed x axis as in long contour plots
source
FLOWFarm.adjust_for_wind_shearMethod
adjust_for_wind_shear(locz, reference_velocity, reference_height, ground_height, model::PowerLawWindShear)

Uses provided velocity at a given height to estimate the velocity at a different height due to wind shear. Ground height may be tuned because the power law does not always hold near the ground.

Arguments

  • locz::Float: height of desired velocity
  • reference_velocity::Float: known velocity at reference_height
  • reference_height::Float: height of known velocity
  • ground_height::Float: height of the ground (typically zero)
  • model::AbstractWindShearModel: wind shear model to use for calculations
source
FLOWFarm.boundary_normals_calculatorMethod
boundary_normals_calculator(boundary_vertices; nboundaries=1)

Outputs the unit vectors perpendicular to each edge of each polygon in a set of polygons, given the Cartesian coordinates for each polygon's vertices.

Arguments

  • boundary_vertices::Array{Float,1} : ragged array of arrays containing all the boundary vertices of each polygon, counterclockwise
  • nboundaries::Int : the number of boundaries in the set
source
FLOWFarm.calcMinorAngleFunction
calcMinorAngle(bndry_x, bndry_y, bndry_z=[0,0,0])

Given three points in space, calculates the magnitude of the non-reflex angle formed at the center point. Created to be used in VRbounarystartup()

Arguments

  • bndry_x::Array{Float,1} : 1-D array of x-coordinates for the vertices around a singlar closed boundary
  • bndry_y::Array{Float,1} : 1-D array of y-coordinates for the vertices around a singlar closed boundary
  • bndry_z::Array{Float,1} : 1-D array of z-coordinates for the vertices around a singlar closed boundary. Default to [0,0,0] for X-Y plane
source
FLOWFarm.calcSmallestAngleMethod
calcSmallestAngle(bndry_x_clsd, bndry_y_clsd)

Given a 1-D closed array of boundary verticies (with first point repeated at the end) it determines the smallest non-reflex angle created by any three consecutive verticies along the boundary. Created to be used in VRbounarystartup()

Arguments

  • bndry_x::Array{Float,1} : 1-D array of x-coordinates for the vertices around a singlar closed boundary
  • bndry_y::Array{Float,1} : 1-D array of y-coordinates for the vertices around a singlar closed boundary
source
FLOWFarm.calc_moment_stressFunction
calc_moment_stress(mx,my,dx,dy,Rcyl=1.771,tcyl=0.06)

Calculates stresses from bending moments on a hollow cylinder

Arguments

  • mx::Float: x moment
  • my::Float: y moment
  • dx::Float: x distance to the location of interest
  • dy::Float: y distance to the location of interest

Keyword Arguments

  • Rcyl::Float: radius of the cylinder
  • tcyl::Float: thickenss of the cylinder
source
FLOWFarm.calculate_aepMethod
calculate_aep(turbine_x, turbine_y, turbine_z, rotor_diameter,
hub_height, turbine_yaw, ct_model, generator_efficiency, cut_in_speed,
cut_out_speed, rated_speed, rated_power, wind_resource, power_models::Array{AbstractPowerModel}, model_set::AbstractModelSet;
rotor_sample_points_y=[0.0], rotor_sample_points_z=[0.0], hours_per_year=365.25*24.0)

Calculate wind farm AEP

Arguments

  • turbine_x::Array{Float,nTurbines}: turbine east-west locations in the global reference frame
  • turbine_y::Array{Float,nTurbines}: turbine north-south locations in the global reference frame
  • turbine_z::Array{Float,nTurbines}: turbine base height in the global reference frame
  • rotor_diameter::Array{Float,nTurbines}
  • hub_height::Array{TF,nTurbines}: turbine hub heights
  • turbine_yaw::Array{TF,nTurbines}: turbine yaw for the given wind direction in radians
  • ct_model::AbstractThrustCoefficientModel: defines how the thrust coefficient changes with state etc
  • generator_efficiency::Array{Float,nTurbines}
  • cut_in_speed::Array{Float,nTurbines}
  • cut_out_speed::Array{Float,nTurbines}
  • rated_speed::Array{Float,nTurbines}
  • rated_power::Array{Float,nTurbines}
  • wind_resource::DiscretizedWindResource: wind resource discreption (directions, speeds, frequencies, etc)
  • power_model::Array{): elements of array should be sub types of AbstractPowerModel
  • model_set::AbstractModelSet: defines wake-realated models to be used in analysis
  • rotor_sample_points_y::Array{TF,N}: horizontal wind location of points to sample across the rotor swept area when calculating the effective wind speed for the wind turbine. Points are centered at the hub (0,0) and scaled by the radius (1=tip of blades)
  • rotor_sample_points_z::Array{TF,N}: vertical wind location of points to sample across the rotor swept area when calculating the effective wind speed for the wind turbine. Points are centered at the hub (0,0) and scaled by the radius (1=tip of blades)
  • hours_per_year::Float: hours per year (averaged for leap year by default)
source
FLOWFarm.calculate_ctMethod
calculate_ct(model::ThrustModelConstantCt)

Calculate the thrust coefficient for a wind turbine based on a pre-determined constant ct

Arguments

  • inflow_velocity::Float: inflow velocity of the wind turbine (unused for const. ct)
  • thrust_model::ThrustModelConstantCt: struct containing a constant ct value for computation
source
FLOWFarm.calculate_ctMethod
calculate_ct(inflow_velocity, thrust_model::ThrustModelCtPoints)

Calculate the thrust coefficient for a wind turbine based on a pre-determined ct curve
    with linear interpolation.

Arguments

  • inflow_velocity::Float: inflow velocity of the wind turbine
  • thrust_model::ThrustModelCtPoints: Struct containing ct and velocity points for ct curve
source
FLOWFarm.calculate_flow_fieldMethod

calculateflowfield(xrange, yrange, zrange, modelset::AbstractModelSet, turbinex, turbiney, turbinez, turbineyaw, turbinect, turbineai, rotordiameter, hubheight, turbinelocalti, sortedturbineindex, wtvelocities, windresource; windfarmstate_id=1)

Generates a flow field for a given state and cross section

Arguments

  • xrange::Range: range defining east-west locations to sample in global reference frame
  • yrange::Range: range defining north-west locations to sample in global reference frame
  • zrange::Range: range defining vertical locations to sample in global reference frame
  • model_set::AbstractModelSet: defines wake-realated models to be used in analysis
  • turbine_x::Array{TF,nTurbines}: turbine east-west locations in the global reference frame
  • turbine_y::Array{TF,nTurbines}: turbine north-south locations in the global reference frame
  • turbine_z::Array{TF,nTurbines}: turbine base height in the global reference frame
  • turbine_yaw::Array{TF,nTurbines}: turbine yaw for the given wind direction in radians
  • turbine_ct::Array{TF,nTurbines}: thrust coefficient of each turbine for the given state
  • turbine_ai::Array{TF,nTurbines}: turbine axial induction for the given state
  • rotor_diameter::Array{TF,nTurbines}: turbine rotor diameters
  • hub_height::Array{TF,nTurbines}: turbine hub heights
  • turbine_local_ti::Array{TF,nTurbines}: turbine local turbulence intensity for the given state
  • sorted_turbine_index::Array{TF,nTurbines}: turbine north-south locations in the global reference frame
  • wtvelocities::Array{TF,nTurbines}: effective inflow wind speed for given state
  • wind_resource::DiscretizedWindResource: wind resource discreption (directions, speeds, frequencies, etc)
  • wind_farm_state_id::Int: index to correct state to use from wind resource provided. Defaults to 1
source
FLOWFarm.calculate_local_tiMethod
calculate_local_ti(turbine_x, turbine_y, ambient_ti, rotor_diameter, hub_height, turbine_yaw, turbine_local_ti, sorted_turbine_index,
turbine_inflow_velcities, turbine_ct, ti_model::LocalTIModelMaxTI; turbine_id=1, tol=1E-6)

Returns local turbulence intensity calculated using Niayifar and Porte Agel 2015, 2016 using smooth max on area TI ratio

Arguments

  • turbine_x::Array{Float,nTurbines}: turbine wind direction locations in the wind direction reference frame
  • turbine_y::Array{Float,nTurbines}: turbine cross wind locations in the wind direction reference frame
  • ambient_ti::Float: ambient turbulence intensity
  • rotor_diameter::Array{Float,nTurbines}: rotor diameters of all turbines
  • hub_height::Array{Float,nTurbines}: hub heights of all turbines relative to the ground
  • turbine_yaw::Array{Float,nTurbines}: yaw of all turbines for the current wind state in radians
  • turbine_local_ti::Array{Float,nTurbines}: local turbulence intensity of all turbines for the current wind state`
  • sorted_turbine_index::Array{Float,nTurbines}: turbine north-south locations in the global reference frame
  • turbine_inflow_velcities::Array{Float,nTurbines}: effective inflow wind speed at each turbine for given state
  • turbine_ct::Array{Float,nTurbines}: thrust coefficient of each turbine for the given state
  • ti_model::LocalTIModelMaxTI: contains a struct defining the desired turbulence intensity model, no local TI in this case
  • turbine_id::Int: index of wind turbine of interest. Provide 1 as default.
  • tol::Float: How far upstream a turbine should be before being included in TI calculations
source
FLOWFarm.calculate_local_tiMethod
calculate_local_ti(turbine_x, turbine_y, ambient_ti, rotor_diameter, hub_height, turbine_yaw, turbine_local_ti, sorted_turbine_index,
turbine_inflow_velcities, turbine_ct, ti_model::LocalTIModelNoLocalTI; turbine_id=1, tol=1E-6)

Returns ambient turbulence intesity value whenever local turbulence intensity is requestesd

Arguments

  • turbine_x::Array{Float,nTurbines}: turbine wind direction locations in the wind direction reference frame
  • turbine_y::Array{Float,nTurbines}: turbine cross wind locations in the wind direction reference frame
  • ambient_ti::Float: ambient turbulence intensity
  • rotor_diameter::Array{Float,nTurbines}: rotor diameters of all turbines
  • hub_height::Array{Float,nTurbines}: hub heights of all turbines relative to the ground
  • turbine_yaw::Array{Float,nTurbines}: yaw of all turbines for the current wind state in radians
  • turbine_local_ti::Array{Float,nTurbines}: local turbulence intensity of all turbines for the current wind state`
  • sorted_turbine_index::Array{Float,nTurbines}: turbine north-south locations in the global reference frame
  • turbine_inflow_velcities::Array{Float,nTurbines}: effective inflow wind speed at each turbine for given state
  • turbine_ct::Array{Float,nTurbines}: thrust coefficient of each turbine for the given state
  • ti_model::LocalTIModelNoLocalTI: contains a struct defining the desired turbulence intensity model, no local TI in this case
  • turbine_id::Int: index of wind turbine of interest. Provide 1 as default.
  • tol::Float: How far upstream a turbine should be before being included in TI calculations
source
FLOWFarm.calculate_powerMethod
calculate_power(generator_efficiency, air_density, rotor_area, wt_velocity, power_model)

Calculate the power for a wind turbine based on standard theory for region 2 using a constant cp

Arguments

  • generator_efficiency::Float: Efficiency of the turbine generator
  • air_density::Float: Air density
  • rotor_area::Float: Rotor-swept area of the wind turbine
  • wt_velocity::Float: Inflow velocity to the wind turbine
  • cut_in_speed::Float: cut in speed of the wind turbine
  • rated_speed::Float: rated speed of the wind turbine
  • cut_out_speed::Float: cut out speed of the wind turbine
  • rated_power::Float: rated power of the wind turbine
  • power_model::PowerModelConstantCp: Struct containing the cp value to be used in region 2
source
FLOWFarm.calculate_powerMethod
calculate_power(generator_efficiency, air_density, rotor_area, wt_velocity, cut_in_speed, rated_speed, cut_out_speed, rated_power, power_model)

Calculate the power for a wind turbine based on a cp curve with linear interpolation

Arguments

  • generator_efficiency::Float: Efficiency of the turbine generator
  • air_density::Float: Air density
  • rotor_area::Float: Rotor-swept area of the wind turbine
  • wt_velocity::Float: Inflow velocity to the wind turbine
  • cut_in_speed::Float: cut in speed of the wind turbine
  • rated_speed::Float: rated speed of the wind turbine
  • cut_out_speed::Float: cut out speed of the wind turbine
  • rated_power::Float: rated power of the wind turbine
  • power_model::PowerModelCpPoints: Struct containing the velocity and cp values defining the cp curve
source
FLOWFarm.calculate_powerMethod
calculate_power(generator_efficiency, air_density, rotor_area, wt_velocity, cut_in_speed, rated_speed, cut_out_speed, rated_power, power_model)

Calculates wind turbine power using a cubic estimation based on turbine specifications as defined in https://github.com/byuflowlab/iea37-wflo-casestudies/blob/master/cs3-4/iea37-cs3-announcement.pdf

Arguments

  • generator_efficiency::Float: Efficiency of the turbine generator
  • air_density::Float: Air density
  • rotor_area::Float: Rotor-swept area of the wind turbine
  • wt_velocity::Float: Inflow velocity to the wind turbine
  • cut_in_speed::Float: cut in speed of the wind turbine
  • rated_speed::Float: rated speed of the wind turbine
  • cut_out_speed::Float: cut out speed of the wind turbine
  • rated_power::Float: rated power of the wind turbine
  • power_model::PowerModelPowerCurveCubic: Empty struct
source
FLOWFarm.calculate_powerMethod
calculate_power(generator_efficiency, air_density, rotor_area, wt_velocity, cut_in_speed, rated_speed, cut_out_speed, rated_power, power_model)

Calculate the power for a wind turbine based on a pre-determined power curve with linear interpolation

Arguments

  • generator_efficiency::Float: Efficiency of the turbine generator
  • air_density::Float: Air density
  • rotor_area::Float: Rotor-swept area of the wind turbine
  • wt_velocity::Float: Inflow velocity to the wind turbine
  • cut_in_speed::Float: cut in speed of the wind turbine
  • rated_speed::Float: rated speed of the wind turbine
  • cut_out_speed::Float: cut out speed of the wind turbine
  • rated_power::Float: rated power of the wind turbine
  • power_model::PowerModelPowerPoints: Struct containing the velocity and power values defining the power curve
source
FLOWFarm.calculate_power_from_cpMethod
calculate_power_from_cp(generator_efficiency, air_density, rotor_area, cp, wt_velocity)

Calculate the power for a wind turbine based on standard theory for region 2

Arguments

  • generator_efficiency::Float: Efficiency of the turbine generator
  • air_density::Float: Air density
  • rotor_area::Float: Rotor-swept area of the wind turbine
  • cp::Float: Power coefficient of the wind turbine
  • wt_velocity::Float: Inflow velocity to the wind turbine
source
FLOWFarm.calculate_state_aepsMethod
calculate_state_aeps(turbine_x, turbine_y, turbine_z, rotor_diameter,
hub_height, turbine_yaw, ct_model, generator_efficiency, cut_in_speed,
cut_out_speed, rated_speed, rated_power, wind_resource, power_models::Array{AbstractPowerModel}, model_set::AbstractModelSet;
rotor_sample_points_y=[0.0], rotor_sample_points_z=[0.0])

Calculate AEP for each requested state respectively

Arguments

  • turbine_x::Array{Float,nTurbines}: turbine east-west locations in the global reference frame
  • turbine_y::Array{Float,nTurbines}: turbine north-south locations in the global reference frame
  • turbine_z::Array{Float,nTurbines}: turbine base height in the global reference frame
  • rotor_diameter::Array{Float,nTurbines}
  • hub_height::Array{TF,nTurbines}: turbine hub heights
  • turbine_yaw::Array{TF,nTurbines}: turbine yaw for the given wind direction in radians
  • ct_model::AbstractThrustCoefficientModel: defines how the thrust coefficient changes with state etc
  • generator_efficiency::Array{Float,nTurbines}
  • cut_in_speed::Array{Float,nTurbines}
  • cut_out_speed::Array{Float,nTurbines}
  • rated_speed::Array{Float,nTurbines}
  • rated_power::Array{Float,nTurbines}
  • wind_resource::DiscretizedWindResource: wind resource discreption (directions, speeds, frequencies, etc)
  • power_model::Array{nTurbines}: elemenst of array should be sub-types of AbstractPowerModel
  • model_set::AbstractModelSet: defines wake-realated models to be used in analysis
  • rotorsamplepoints_y::Array{TF,N}`: horizontal wind location of points to sample across the rotor swept area when calculating the effective wind speed for the wind turbine. Points are centered at the hub (0,0) and scaled by the radius (1=tip of blades)
  • rotorsamplepoints_z::Array{TF,N}`: vertical wind location of points to sample across the rotor swept area when calculating the effective wind speed for the wind turbine. Points are centered at the hub (0,0) and scaled by the radius (1=tip of blades)
  • hours_per_year::Float: hours per year (averaged for leap year by default)
source
FLOWFarm.calculate_state_turbine_powersMethod
calculate_state_turbine_powers(turbine_x, turbine_y, turbine_z, rotor_diameter,
hub_height, turbine_yaw, ct_model, generator_efficiency, cut_in_speed,
cut_out_speed, rated_speed, rated_power, wind_resource, power_models, model_set::AbstractModelSet;
rotor_sample_points_y=[0.0], rotor_sample_points_z=[0.0])

Calculate power for each turbine for all states respectively

Arguments

  • turbine_x::Array{Float,nTurbines}: turbine east-west locations in the global reference frame
  • turbine_y::Array{Float,nTurbines}: turbine north-south locations in the global reference frame
  • turbine_z::Array{Float,nTurbines}: turbine base height in the global reference frame
  • rotor_diameter::Array{Float,nTurbines}
  • hub_height::Array{TF,nTurbines}: turbine hub heights
  • turbine_yaw::Array{TF,nTurbines}: turbine yaw for the given wind direction in radians
  • ct_model::AbstractThrustCoefficientModel: defines how the thrust coefficient changes with state etc
  • generator_efficiency::Array{Float,nTurbines}
  • cut_in_speed::Array{Float,nTurbines}
  • cut_out_speed::Array{Float,nTurbines}
  • rated_speed::Array{Float,nTurbines}
  • rated_power::Array{Float,nTurbines}
  • wind_resource::DiscretizedWindResource: wind resource discreption (directions, speeds, frequencies, etc)
  • power_models::Array{nTurbines}: elemenst of array should be sub-types of AbstractPowerModel
  • model_set::AbstractModelSet: defines wake-realated models to be used in analysis
  • rotorsamplepoints_y::Array{TF,N}`: horizontal wind location of points to sample across the rotor swept area when calculating the effective wind speed for the wind turbine. Points are centered at the hub (0,0) and scaled by the radius (1=tip of blades)
  • rotorsamplepoints_z::Array{TF,N}`: vertical wind location of points to sample across the rotor swept area when calculating the effective wind speed for the wind turbine. Points are centered at the hub (0,0) and scaled by the radius (1=tip of blades)
source
FLOWFarm.calculate_turbine_powerMethod
calculate_turbine_power(generator_efficiency, cut_in_speed, cut_out_speed, rated_speed, 
rated_power, rotor_diameter, wt_velocity, power_model::AbstractPowerModel, air_density)

Calculate the power for all wind turbines. Dispaches to desired power model.

Arguments

  • generator_efficiency::Array{Float,nTurbines}
  • cut_in_speed::Array{Float,nTurbines}
  • cut_out_speed::Array{Float,nTurbines}
  • rated_speed::Array{Float,nTurbines}
  • rated_power::Array{Float,nTurbines}
  • rotor_diameter::Array{Float,nTurbines}
  • wt_velocity::Array{Float,nTurbines}: turbine effective wind speeds for current state only
  • `power_model::AbstractPowerModel)
  • air_density::Float
source
FLOWFarm.circle_boundaryMethod
circle_boundary(center,radius,turbine_x,turbine_y)

calculate the distance from each turbine to a circular boundary. Negative means the turbine is inside the boundary

Arguments

  • center::Float: circular boundary center [x,y]
  • radius::Float: circulat boundary radius
  • turbine_x::Array{Float}: turbine x locations
  • turbine_y::Array{Float}: turbine y locations
source
FLOWFarm.closeBndryListMethod
closeBndryList(bndryPts_x, bndryPts_y)

Appends the 1st element to the end of the coordinate arrays if it is not already repeated. Note, this will only work on 1-D arrays. For an array of 1-D arrays, use closeBndryLists(bndryPts_x, bndryPts_y) (note the plural, not singular 'Lists' in the function title)

Arguments

  • bndryPts_x::Array{Float,1} : 1-D array of x-coordinates for the vertices around a singlar closed boundary
  • bndryPts_y::Array{Float,1} : 1-D array of y-coordinates for the vertices around a singlar closed boundary
source
FLOWFarm.closeBndryListsMethod
closeBndryLists(bndryPts_x, bndryPts_y)

Appends the 1st element to the end of each array for a closed boundary. Note, this will not function properly if there is only one region. For only one region, use closeBndryList(bndryPts_x, bndryPts_y) (note the singular, not plural 'List' in the function title)

Arguments

  • bndryPts_x::Array{Float,1} : N-D array of x-coordinates for the vertices around N-many closed boundaries
  • bndryPts_y::Array{Float,1} : N-D array of y-coordinates for the vertices around N-many closed boundaries
source
FLOWFarm.convex_boundaryMethod
convex_boundary(boundary_vertices,boundary_normals,turbine_x,turbine_y)

calculate the distance from each turbine to a possibly non-circular, but convex boundary. Negative means the turbine is inside the boundary

Arguments

  • boundary_vertices::Array{Float,2}: vertices of the convex hull CCW in order s.t. boundaryVertices[i] -> first point of face for unit_normals[i]
  • boundary_normals::Array{Float,2}: unit normal vector for each boundary face CCW where boundaryVertices[i] is the first point of the corresponding face
  • turbine_x::Array{Float}: turbine x locations
  • turbine_y::Array{Float}: turbine y locations
source
FLOWFarm.coordDistMethod
coordDist(x1, y1, x2, y2)

Given a two points (x1, y1) and (x2, y2), returns the euclidean distance between them

Arguments

  • x1::Float64 : x-coord of the first point
  • y1::Float64 : y-coord of the first point
  • x2::Float64 : x-coord of the second point
  • y2::Float64 : y-coord of the second point
source
FLOWFarm.distributed_velocity_opFunction
distributed_velocity_op(V, Omega, r, precone, yaw, tilt, azimuth, rho, mu=1.81206e-05, asound=1.0)

Return the operating points along the blade considering varied inflow along the blade.

Arguments

  • V::Array{Float}: velocity inflow at each r
  • Omega::Float: rotor rotational speed
  • r::Array{Float}: radial locations of interest
  • precone::Float: rotor precone angle
  • yaw::Float: rotor yaw angle
  • tilt::Float: rotor tilt angle
  • azimuth::Float: blade azimuth angle
  • rho::Float: air density

Keyword Arguments

  • mu::Float: air viscocity (can usually use the default)
  • asound::Float: speed of sound (can usually use the default)
source
FLOWFarm.find_upstream_turbinesMethod
find_upstream_turbines(turbinex, turbiney, winddirection, diameter; inverse=false)

A convenience function to quickly find either which turbines are waked, or those that are not.

Arguments

  • turbinex::Array{T,1}: x locations of turbines in global reference frame
  • turbiney::Array{T,1}: y locations of turbines in global reference frame
  • winddirection::Real or winddirection::AbstractArray: wind direction in radians in meteorological coordinates (0 rad. = from North)
  • diameter::Array{T,1}: diameters of all wind turbines
source
FLOWFarm.find_xyz_simpleMethod
find_xyz_simple(x_hub,y_hub,z_hub,r,yaw,azimuth)

Find the xyz locations of points along a blade given it's location and azimuth angle. Currently doesn't consider precone or tilt.

Arguments

  • x_hub::Float: x location of hub
  • y_hub::Float: y location of hub
  • z_hub::Float: z location of hub (hub height if no topology)
  • r::Array{Float}: radial locations of interest
  • precone::Float: rotor precone angle
  • yaw::Float: rotor yaw angle
  • azimuth::Float: blade azimuth angle
source
FLOWFarm.getNextFileNameFunction
getNextFileName(directory, file_name, file_type)

Checks if a file of the given directory and name exists. If not, increments to the next index so as not to overwrite previously written files. If it reaches the max number of overwrites, it will default to <directory/filename.filetype> To default to this, set <max_check=0> in function call.

Arguments

  • directory::String: path/to/write/file/at/
  • file_name::String: Whatever the name of the file desired
  • file_type::String: ex "yaml", "txt", "csv", etc...
  • max_check::Int: the maximum number of files to check
source
FLOWFarm.getPerimeterLengthMethod
getPerimeterLength(bndry_x_clsd, bndry_y_clsd)

Given a 1-D closed array of boundary verticies (with first point repeated at the end) returns the length along the perimeter. Created to be used in VRbounarystartup()

Arguments

  • bndry_x::Array{Float,1} : 1-D array of x-coordinates for the vertices around a singlar closed boundary
  • bndry_y::Array{Float,1} : 1-D array of y-coordinates for the vertices around a singlar closed boundary
source
FLOWFarm.getUpDwnYvalsMethod
getUpDwnYvals(turbine_x, bndry_x_clsd, bndry_y_clsd, bndry_corner_indcies)

Supplements FLOWFarm's splinedboundary() function by calculating (for a given x location) the maximum and minimum y-value permitted to remain "inside" the boundary. If turbinex is located left of the boundary's leftmost vertex or right of the boundary's rightmost vertex, it return's that corresponding vertex's y-value as the max and min, as default. Returns two values, the minimum and maximum interior y-values withing a boundary for the given turbine_x value. Note that all boundary coordinates must be in the first quadrant of the Cartesian coordinate system (+x and +y values only)

Arguments

  • turbine_x::Array{Float}: x-value of the turbine being examined
  • bndry_x_clsd::Array{Float}: x locations of boundary vertices, CCW in order s.t. boundaryVertices[0] is the NE corner of boundary, and with the first vertex duplicated at the end for completeness of calcs.
  • bndry_y_clsd::Array{Float}: y locations of boundary vertices, CCW in order s.t. boundaryVertices[0] is the NE corner of boundary, and with the first vertex duplicated at the end for completeness of calcs.
  • bndry_corner_indcies::Array{Float}: An array of 3 or 4 indicies in the bndryx/yclsd arrays that correspond to the three four "corners" used between splined "sides"
source
FLOWFarm.get_boundary_yamlMethod
get_boundary_yaml(filename)

Returns the boundaries of a wind farm as defined in a yaml file in the format used in FLOWFarm. Returns N by 2 array for single region farm and an array of N by 2 arrays for multiple regions. Returned regions are sorted alphabetically by the keys provided in the yaml file.

Arguments

  • file_name::String: relative/path/to/file
source
FLOWFarm.get_momentsMethod
get_moments(out,Rhub,Rtip,r,az,precone,tilt)

Trapezoidal integration to find the blade root bending moment using the loads distribution

Arguments

  • `out::CCBlade dict: output from running CCBlade solve
  • Rhub::Float: radius of the rotor hub
  • Rtip::Float: radius of the blade tip
  • r::Array{Float}: radial locations of interest
  • az::Float: blade azimuth angle
  • precone::Float: rotor precone angle
  • tilt::Float: rotor tilt angle
source
FLOWFarm.get_peaksMethod
get_peaks(array)

get the turning point values of a signal

Arguments

  • array::Array{Float}: the signal to find the turning points, or peaks
source
FLOWFarm.get_peaks_indicesMethod
get_peaks_indices(array)

return the indices of the signal peaks

Arguments

  • array::Array{Float}: the signal to find the turning points, or peaks
source
FLOWFarm.get_turb_loc_YAMLMethod
get_turb_loc_YAML(file_name)

read in turbine locations and related problem file names from .yaml

Arguments

  • file_name::String: path/and/name/of/location/file.yaml
source
FLOWFarm.get_wind_rose_YAMLMethod
get_wind_rose_YAML(file_name)

read in wind resource information from .yaml

Arguments

  • file_name::String: path/to/wind/resource/file.yaml
source
FLOWFarm.grid_pointsMethod
grid_points(n)

Generates points in a grid. If n is not a perfect square, then the nearest square root will be used for the side length of the grid.

Arguments

  • n::Float: number of points to generate
source
FLOWFarm.hermite_splineMethod
hermite_spline(x, x0, x1, y0, dy0, y1, dy1)

Produces the y and (optionally) dy values for a hermite cubic spline interpolating between two end points with known slopes

Arguments

  • x::Float: x position of output y
  • x0::Float: x position of upwind endpoint of spline
  • x1::Float: x position of downwind endpoint of spline
  • y0::Float: y position of upwind endpoint of spline
  • dy0::Float: slope at upwind endpoint of spline
  • y1::Float: y position of downwind endpoint of spline
  • dy1::Float: slope at downwind endpoint of spline
source
FLOWFarm.iea37cs4BndryVRIntPMMethod
iea37cs4BndryVRIntPM(bndry_x_clsd, bndry_y_clsd, bndry_corner_indcies, turbine_x, turbine_y, turb_diam, turb_min_space, num_turbs_to_place)

Uses the Variable reduction method for placing boundary turbines, and the Partition Method (from splined_boundary()) for random interior points, maintaining proper spacing from all previously placed turbines.

Arguments

  • bndry_x_clsd::Array{Float,1} : 1-D array of x-coordinates for the vertices around a singlar closed boundary
  • bndry_y_clsd::Array{Float,1} : 1-D array of y-coordinates for the vertices around a singlar closed boundary
  • bndry_corner_indcies::Float64: The indicies within <bndryxclsd> and <bndryyclsd> which denote the "corners" adjacent turbines
  • 'turbminspace::Float64`: For proximity knowledge, the minimum spacing required between any two turbines
  • 'numbndryturbs::Float64`: The number of turbines desired to be placed along the boundary. If too many are selected (due to spacing condtraints), the remaining will be placed in the interior
  • 'numtotturbs::Float64`: The number of total turbines to be placed both on the boundary and in the interior
source
FLOWFarm.latlong_to_xyMethod
latlong_to_xy(latitude, longitude, utm_zone; isnorth=true, units="m")

Converts arrays of points from latitude and longitude to x and y in meters in a local coordinate frame based on the point with the lowest magnitude latitude,

Arguments

  • latitude::Array{Float,N}
  • longitude::Array{Float,N}
  • isnorth::Float: specifies if the point is in the northern hemisphere (defaul: true)
source
FLOWFarm.met2cartMethod
met2cart(angle_met)

Convert from meteorological polar system (CW, 0 rad.=N, wind from) to cartesian polar system (CCW, 0 rad.=E, wind to).

Arguments

  • angle_met::Number: an angle in radians in a meteorological coordinate system
source
FLOWFarm.multiple_components_opFunction
multiple_components_op(U, V, W, Omega, r, precone, yaw, tilt, azimuth, rho, mu=1.81206e-05, asound=1.0)

Return the operating points along the blade considering all the inflow velocity components.

Arguments

  • U::Array{Float}: u velocity component of the inflow at each r
  • V::Array{Float}: v velocity component of the inflow at each r
  • W::Array{Float}: w velocity component of the inflow at each r
  • Omega::Float: rotor rotational speed
  • r::Array{Float}: radial locations of interest
  • precone::Float: rotor precone angle
  • yaw::Float: rotor yaw angle
  • tilt::Float: rotor tilt angle
  • azimuth::Float: blade azimuth angle
  • rho::Float: air density

Keyword Arguments

  • mu::Float: air viscocity (can usually use the default)
  • asound::Float: speed of sound (can usually use the default)
source
FLOWFarm.nansafenormMethod
nansafenorm(v)

Calculate the norm of a vector, but if the sum of the squares is less than the given tolerance then use the line y = a(sqrt(eps())/eps()) so that the derivative is well defined.

Arguments

  • v::Vector{}: takes the norm of this vector, but avoids NaN by using a linear approximation of sqrt near 0.
source
FLOWFarm.nansafesqrtMethod
nansafesqrt(a)

Calculate the square root of a number, but if the number is less than the given tolerance then use the line y = a(sqrt(eps())/eps()) so that the derivative is well defined.

Arguments

  • a::Number: takes the square root of this value, or approximates it with a line for a < eps()
source
FLOWFarm.plotboundary!Method
plotboundary!(ax, boundary_vertices; color="k", linestyle="-", nboundaries=1)

Convenience function for plotting wind farm boundaries

Arguments

  • `ax
  • boundary_vertices::Array{Float,1}(nvertices): an nx2 array of boundary vertices for polygon or [[centerx, centery], radius] for circle boundary
  • color::=String: sets color for turbine markers
  • linestyle::String: sets the line style for the boundary
  • nboundaries::Int: number of discrete boundary regions
source
FLOWFarm.plotlayout!Method
plotlayout!(ax, turbinex, turbiney, rotordiameter; aspect="equal", xlim=[], ylim=[], fill=false, color="k", markeralpha=1, title="")

Convenience function for plotting wind farm layouts

Arguments

  • `ax
  • turbinex::Array{Float,1}(nturbines): an array x coordinates of wind turbine locations
  • turbiney::Array{Float,1}(nturbines): an array y coordinates of wind turbine locations
  • rotordiameter::Array{Float,1}(nturbines): an array rotor diameters of wind turbines
  • aspect::String: set plot aspect ratio, default="equal"
  • xlim::Array: limits in x coordinate. "[]" results in limits being automatically defined
  • ylim::Array: limits in y coordinate. "[]" results in limits being automatically defined
  • fill::Bool: determines whether turbine circle markers are filled or not
  • color::=String: sets color for turbine markers
  • markeralpha::Int: determines tranparancy of turbine markers
  • itle::String: optional title to include on the plot
source
FLOWFarm.plotrotorsamplepoints!Method
plotrotorsamplepoints!(ax, y, z; rotordiameter=2.0, aspect="equal", ylim=[], zlim=[], fill=false, color="k", markeralpha=1, title="")

Convenience function for plotting where points are being sampled on the wind turbine rotor

Arguments

  • `ax
  • y::Array{Float,1}(nturbines): an array x coordinates of wind turbine locations
  • z::Array{Float,1}(nturbines): an array y coordinates of wind turbine locations
  • rotordiameter::Number: rotor diameter of wind turbine
  • aspect::String: set plot aspect ratio, default="equal"
  • ylim::Array: limits in y coordinate. "[]" results in limits being automatically defined
  • zlim::Array: limits in z coordinate. "[]" results in limits being automatically defined
  • fill::Bool: determines whether turbine circle markers are filled or not
  • color::=String: sets color for turbine markers
  • markeralpha::Int: determines tranparancy of turbine markers between 0 (transparent) and 1 (opaque)
  • itle::String: optional title to include on the plot
source
FLOWFarm.plotsingleboundary!Method
plotsingleboundary!(ax, boundary_vertices; color="k", linestyle="-')

Convenience function for plotting wind farm boundaries

Arguments

  • `ax
  • boundary_vertices::Array{Float,1}(nvertices): an nx2 array of boundary vertices for polygon or [[centerx, centery], radius] for circle boundary
  • color::=String: sets color for turbine markers
  • linestyle::String: sets the line style for the boundary
source
FLOWFarm.plotwindfarm!Method
plotwindfarm!(ax, boundary_vertices, turbinex, turbiney, rotordiameter; nboundaries=1, aspect="equal", xlim=[], ylim=[], fill=false, color="k", markeralpha=1, title="")

Convenience function for plotting wind farms

Arguments

  • `ax
  • boundary_vertices::Array{Float,1}(nvertices): an nx2 array of boundary vertices
  • turbinex::Array{Float,1}(nturbines): an array x coordinates of wind turbine locations
  • turbiney::Array{Float,1}(nturbines): an array y coordinates of wind turbine locations
  • rotordiameter::Array{Float,1}(nturbines): an array rotor diameters of wind turbines
  • nboundaries::Int: number of discrete boundary regions
  • aspect::String: set plot aspect ratio, default="equal"
  • xlim::Array: limits in x coordinate. "[]" results in limits being automatically defined
  • ylim::Array: limits in y coordinate. "[]" results in limits being automatically defined
  • fill::Bool: determines whether turbine circle markers are filled or not
  • color::=String: sets color for turbine markers
  • markeralpha::Int: determines tranparancy of turbine markers
  • title::String: optional title to include on the plot
source
FLOWFarm.plotwindresource!Method
plotwindresource!(ax::Array, windresource::ff.DiscretizedWindResource; roundingdigits=[1,3], fill=false, alpha=0.5, colors=["b", "b"], fontsize=8, edgecolor=nothing, rlabel_position=-45)

Convenience function for visualizing the wind speed and wind frequency roses

Arguments

  • ax::Array: pre-initialized single dimension array of axes from pyplot with length at least 2
  • windresource::ff.DiscretizedWindResource: wind rose information
  • roundingdigits::Array{Int, 1}: how many significant digits to round to on each axis in ax
  • fill::Bool: determines whether bars are filled or not
  • alpha::Number: tranparancy of bars between 0 (transparent) and 1 (opaque)
  • colors::=Array{String, 1}: sets color for turbine markers
  • fontsize::Int: font size of text on figures
  • edgecolor: color of edges of each bar in polar chart, nothing means no color
  • rlabel_position:Number: Angle at which to draw the radial axes
source
FLOWFarm.plotwindrose!Method
plotwindrose!(ax, d, f; roundingdigit=1, color="C0",alpha=0.5,fontsize=8,
dticks=(0,pi/4,pi/2,3*pi/4,pi,5*pi/4,3*pi/2,7*pi/4),
dlabels=("E","NE","N","NW","W","SW","S","SE"),
fticks=nothing, flabels=nothing, normalize=false, edgecolor=nothing, units="",
rlabel_position=-45)

Convenience function for creating a windrose from any polar data

Arguments

  • ax::PyCall.PyObject: pre-initialized axis from pyplot
  • d::Vector: wind rose directions
  • f::Vector: wind rose radial variable
  • roundingdigit::Int: how many significant digits to round to
  • fontsize::Int: font size of text on figures
  • dticks::Tuple: contains angular tick locations
  • dlabels::Tuple: contains angular tick labels
  • fticks::Tuple: contains radial tick locations
  • flabels::Tuple: contains radial tick labels
  • normalize::Bool: choose whether or not to normalize by the sum of f
  • units::String: Units to append to flabels
  • rlabel_position:Number: Angle at which to draw the radial axis
  • plotcommand::String: Type of plot. Can be ["bar", "plot"]
  • kwargs::Tuple: tuple containing key word arguments to plotcommand in the form (key => "value")
source
FLOWFarm.point_velocityMethod
point_velocity(loc, turbine_x, turbine_y, turbine_z, turbine_yaw, turbine_ct, turbine_ai,
rotor_diameter, hub_height, turbine_local_ti, sorted_turbine_index, wtvelocities,
wind_resource, model_set::AbstractModelSet;
wind_farm_state_id=1, downwind_turbine_id=0)

Calculates the wind speed at a given point for a given state

Arguments

  • loc::Array{TF,3}: Location of interest
  • turbine_x::Array{TF,nTurbines}: turbine east-west locations in the state reference frame
  • turbine_y::Array{TF,nTurbines}: turbine north-south locations in the state reference frame
  • turbine_z::Array{TF,nTurbines}: turbine base height in the state reference frame
  • turbine_yaw::Array{TF,nTurbines}: turbine yaw for the given wind direction in radians
  • turbine_ct::Array{TF,nTurbines}: turbine thrust coefficients for the given state
  • turbine_ai::Array{TF,nTurbines}: turbine axial induction for the given state
  • rotor_diameter::Array{TF,nTurbines}: turbine rotor diameters
  • hub_height::Array{TF,nTurbines}: turbine hub heights
  • turbine_local_ti::Array{TF,nTurbines}: turbine local turbulence intensity for the given state
  • sorted_turbine_index::Array{TF,nTurbines}: array containing indices of wind turbines from most upwind to most downwind turbine in the given state
  • wtvelocities::Array{TF,nTurbines}: effective inflow wind speed for given state
  • wind_resource::DiscretizedWindResource: contains wind resource discreption (directions, speeds, frequencies, etc)
  • wind_farm_state_id::Int: index to correct state to use from wind resource provided. Defaults to 1
  • downwind_turbine_id::Int: index of wind turbine of interest (if any). If not a point for calculating effective wind speed of a turbine, then provide 0 (default)
source
FLOWFarm.pointinpolygonFunction
pointinpolygon(point, vertices, normals=nothing; s=700, method="raycasting", shift=1E-10, return_distance=true)

Given a polygon determined by a set of vertices, determine the signed distance from the point to the polygon.

Returns the negative (-) distance if the point is inside or on the polygon, positive (+) otherwise. If return_distance is set to false, then returns -1 if in polygon or on the boundary, and 1 otherwise.

Arguments

  • point::Vector{Number}(2): point of interest
  • vertices::Vector{Matrix{Number}(2): vertices of polygon
  • normals::Vector{Matrix{Number}(2): if not provided, they will be calculated
  • s::Number: smoothing factor for ksmax function (smoothmax)
  • method::String: currently only raycasting is available
  • shift::Float: how far to shift point if it lies on an edge or vertex
  • return_distance::Bool: if true, return distance. if false, return -1 if in polygon or on the boundary, and 1 otherwise.
source
FLOWFarm.pointonlineMethod
pointonline(p, v1, v2; tol=1E-6)

Given a line determined two points (v1 and v2) determine if the point (p) lies on the line between those points.

Returns true if the point lies on the line (within the given tolerance), false otherwise.

Arguments

  • p::Vector{Number}(2): point of interest
  • v1::Vector{Number}(2): first vertex of the line
  • v2::Vector{Number}(2): second vertex of the line
  • tol::Number: how close the cumulative distance from v1 to p to v2 must be to the distance from v1 to v2 to count as being co-linear
source
FLOWFarm.print_layout_in_cartesian_frame_excelMethod
print_state_layouts_in_cartesian_frame(turbinex, turbiney, winddirections)

Given a wind farm layout in the global reference frame, print the layout rotated to the cartesian frame with wind to the positive x axis (right) for all wind directions.

Arguments

  • turbinex::Array{T,1}: x locations of turbines in global reference frame
  • turbiney::Array{T,1}: y locations of turbines in global reference frame
  • winddirections::Array{T,1}: all wind directions in radians in meteorological coordinates (0 rad. = from North)
source
FLOWFarm.rainflowFunction
rainflow(array_ext,uc_mult=0.5)

Rainflow counting of a signal's turning points

Arguments

    array_ext (numpy.ndarray): array of turning points

Keyword Arguments

    uc_mult (float): partial-load scaling [opt, default=0.5]

Returns

    array_out (numpy.ndarray): (3 x n_cycle) array of rainflow values:
                                1) load range
                                2) range mean
                                3) cycle count
source
FLOWFarm.ray_casting_boundaryMethod
ray_casting_boundary(boundary_vertices,boundary_normals,turbine_x,turbine_y)

Calculate the distance from each turbine to the nearest point on the boundary using the ray-casting algorithm. Negative means the turbine is inside the boundary.

Arguments

  • boundary_vertices::Array{Float,2}: vertices of the boundary CCW in order s.t. boundaryVertices[i] -> first point of face for unit_normals[i]
  • boundary_normals::Array{Float,2}: unit normal vector for each boundary face CCW where boundaryVertices[i] is the first point of the corresponding face
  • turbine_x::Array{Float}: turbine x locations
  • turbine_y::Array{Float}: turbine y locations
  • discrete::Bool: if true, indicates the boundary is made of multiple discrete regions
  • s::Number: smoothing factor for smooth max (ksmax)
  • tol::Float: how close points have to be to vertex or face before they will be shifted slightly to avoid a discontinuity
  • return_region::bool: if true, return a vector specifying which region each turbine is in
source
FLOWFarm.rediscretize_windroseMethod
rediscretize_windrose(windrosein::DiscretizedWindResource, ndirectionbins, nspeedbins)

Function for re-interpreting a wind rose into a desired number of bins. Returns the new
wind rose. Currently only works for windroses with a single speed in each direction.

Arguments

  • windrosein::DiscretizedWindResource: original wind rose
  • ndirectionbins::Integer: number of direction bins for the new wind rose
  • start::Float: direction for first bin in radians
  • averagespeed::Bool: set whether or not to return the average wind speed as the speed for all bins
source
FLOWFarm.rotate_to_wind_directionMethod
rotate_to_wind_direction(xlocs, ylocs, wind_direction_met)

Rotates wind farm coordinates to be in wind direction reference where wind direction is to the positive x.

Arguments

  • xlocs::Array: contains turbine east-west locations in the global reference frame
  • ylocs::Array: contains turbine north-south locations in the global reference frame
  • wind_direction_met::Array: contains wind direction in radians in meteorological standard system (N=0 rad, proceeds CW, wind from direction given)
source
FLOWFarm.rotor_sample_pointsFunction
rotor_sample_points(nsamplepoints=1)

Initializes the sampling locations in the rotor-swept-area. Returns values such that zero is at the turbine hub location and 1 is at the tip of the blades. If a single sample is requested, it will be at the hub location. Otherwise, the points will be located using the sunflower packcing algorithm.

Arguments

  • nsamplepoints::Int: controlls how many sample points to generate
  • alpha::Float: Controls smoothness of the sunflower algorithm boundary. alpha=0 is the standard "jagged edge" sunflower algoirthm and alpha=1 results in a smooth boundary.
  • pradius::Float: the percent of the rotor radius to use in generating initial point grid
  • use_perimeter_points: whether or not to include point exactly on the perimeter of the rotor swept area
source
FLOWFarm.round_farm_random_startMethod
round_farm_random_start(rotor_diameter, center, radius; min_spacing=2., min_spacing_random=3., method="individual")

Generates starting locations for multi-start optimization approaches when the farm boundary is round.

Arguments

  • rotor_diameter::Number: wind turbine diameter
  • center::Number: wind farm center
  • radius::Number: wind farm radius
  • diameter::Array{T,1}: diameters of all wind turbines
source
FLOWFarm.smooth_maxMethod
smooth_max_ndim(x; s=100.0)

Calculate the smooth-max (a.k.a. softmax or LogSumExponential) of the elements in x.

Based on John D. Cook's writings at (1) https://www.johndcook.com/blog/2010/01/13/soft-maximum/ and (2) https://www.johndcook.com/blog/2010/01/20/how-to-compute-the-soft-maximum/

Arguments

  • x::Float: first value for comparison
  • y::Float: second value for comparison
  • s::Float : controls the level of smoothing used in the smooth max
source
FLOWFarm.smooth_maxMethod
smooth_max(x; s=10.0)

Calculate the smoothmax (a.k.a. softmax or LogSumExponential) of the elements in x.

Based on John D. Cook's writings at (1) https://www.johndcook.com/blog/2010/01/13/soft-maximum/ and (2) https://www.johndcook.com/blog/2010/01/20/how-to-compute-the-soft-maximum/

And based on article in FeedlyBlog (3) https://blog.feedly.com/tricks-of-the-trade-logsumexp/

Arguments

  • x::Array{Float,1} : vector with all the input values
  • s::Float : controls the level of smoothing used in the smooth max
source
FLOWFarm.splined_boundaryMethod
splined_boundary(turbine_x, turbine_y, bndry_x_clsd, bndry_y_clsd, bndry_corner_indcies)

calculate the distance from each turbine to a closed boundary made up of zero or more reflex angles (concavities). Boundary will have three or four user-selected "corners", such that the "sides" between corners (that will be splined) are injective functions (meaning that for every x-coord, there exists only one corresponding y-coord). Returns four values for every turbine, corresponding to the distance from the turb to the upper, lower, left, and right splined "sides". A negative return value means the turb is inside the boundary for that "side". Returns a single array of {Float64} of length {length(turbine_x) * 4}. Note that all boundary coordinates must be in the first quadrant of the Cartesian coordinate system (+x and +y values only)

Arguments

  • turbine_x::Array{Float}: turbine x locations
  • turbine_y::Array{Float}: turbine y locations
  • bndry_x_clsd::Array{Float}: x locations of boundary vertices, CCW in order s.t. boundaryVertices[0] is the NE corner of boundary, and with the first vertex duplicated at the end for completeness of calcs.
  • bndry_y_clsd::Array{Float}: y locations of boundary vertices, CCW in order s.t. boundaryVertices[0] is the NE corner of boundary, and with the first vertex duplicated at the end for completeness of calcs.
  • bndry_corner_indcies::Array{Float}: An array of 3 or 4 indicies in the bndryx/yclsd arrays that correspond to the three four "corners" used between splined "sides"
source
FLOWFarm.splined_boundary_discreet_regionsMethod
splined_boundary_discreet_regions(turbine_x, turbine_y, bndry_x_clsd, bndry_y_clsd, bndry_corner_indcies, turbs_per_region)

Uses FLOWFarm's splinedboundary() function to calculate the turbine-boundary constraints for one or more discreet regions, with pre-allocated turbines for each region. Returns four values for every turbine, corresponding to the distance from each turb to the upper, lower, left, and right splined "sides" for the region to which it was allocated. A negative return value means the turb is outside the "side" of boundary for which it has been allocated. Returns a single array of {Float64} of length {length(turbinex) * 4}. Note that all boundary coordinates must be in the first quadrant of the Cartesian coordinate system (+x and +y values only)

Arguments

  • turbine_x::Array{Float}: turbine x locations
  • turbine_y::Array{Float}: turbine y locations
  • bndry_x_clsd::Array{Float}: x locations of boundary vertices, CCW in order s.t. boundaryVertices[0] is the NE corner of boundary, and with the first vertex duplicated at the end for completeness of calcs.
  • bndry_y_clsd::Array{Float}: y locations of boundary vertices, CCW in order s.t. boundaryVertices[0] is the NE corner of boundary, and with the first vertex duplicated at the end for completeness of calcs.
  • bndry_corner_indcies::Array{Float}: An array of 3 or 4 indicies in the bndryx/yclsd arrays that correspond to the three four "corners" used between splined "sides"
  • 'turbsperregion::Array{Int}`: An array of length equivalent to the number of discrete boundary regions, with each element denoting howmany turbines are apportioned to the corresponding region. sum(turbsperregion) must be equivalent to the total number of turbines in the windfarm
source
FLOWFarm.sunflower_pointsMethod
sunflower_points(n; alpha=0.0)

Generates points in a circle of radius=1 using the sunflower packing algorithm.

Arguments

  • n::Float: number of points to generate
  • alpha::Float: Controls the smoothness of the boundary. alpha=0 is the standard "jagged edge" sunflower algoirthm and alpha=1 results in a smooth boundary.
source
FLOWFarm.turbine_powers_one_directionMethod
turbine_powers_one_direction((generator_efficiency, cut_in_speed, cut_out_speed, 
    rated_speed, rated_power, rotor_diameter, turbine_inflow_velcities, air_density, power_model::AbstractPowerModel)

Calculate the power for all wind turbines for a given state

Arguments

  • generator_efficiency::Array{Float,nTurbines}
  • cut_in_speed::Array{Float,nTurbines}
  • cut_out_speed::Array{Float,nTurbines}
  • rated_speed::Array{Float,nTurbines}
  • rated_power::Array{Float,nTurbines}
  • rotor_diameter::Array{Float,nTurbines}
  • turbine_inflow_velcities::Array{Float,nTurbines}: for current state only
  • air_density::Float
  • power_models::Array{nturbines}) elements of array should be be of sub-types or AbstractPowerModel
source
FLOWFarm.turbine_spacingMethod
turbine_spacing(turbine_x,turbine_y)

Calculate the distance between turbines in a wind farm. There is an infinite gradient of this function if two points are exactly the same. This can be avoided by returning the square of the turbine spacing rather than the actual distance, but it makes the gradients scale much more poorly. Because it is very vanishingly rare to have turbines exactly in the same location, this function leaves the square root in the calculations.

Arguments

  • turbine_x::Array{Float}: turbine x locations
  • turbine_y::Array{Float}: turbine y locations
source
FLOWFarm.turbine_velocities_one_directionMethod
point_velocity(turbine_x, turbine_y, turbine_z, rotor_diameter, hub_height, turbine_yaw,
sorted_turbine_index, ct_model, rotor_sample_points_y, rotor_sample_points_z, wind_resource,
model_set::AbstractModelSet; wind_farm_state_id=1)

Calculates the wind speed at a given point for a given state

Arguments

  • turbine_x::Array{TF,nTurbines}: turbine east-west locations in the state reference frame
  • turbine_y::Array{TF,nTurbines}: turbine north-south locations in the state reference frame
  • turbine_z::Array{TF,nTurbines}: turbine base height in the state reference frame
  • rotor_diameter::Array{TF,nTurbines}: turbine rotor diameters
  • hub_height::Array{TF,nTurbines}: turbine hub heights
  • turbine_yaw::Array{TF,nTurbines}: turbine yaw for the given wind direction in radians
  • sorted_turbine_index::Array{TF,nTurbines}: turbine sorted order upstream to downstream for given state
  • ct_model::AbstractThrustCoefficientModel: defines how the thrust coefficient changes with state etc
  • rotorsamplepoints_y::Array{TF,N}`: horizontal wind location of points to sample across the rotor swept area when calculating the effective wind speed for the wind turbine. Points are centered at the hub (0,0) and scaled by the radius (1=tip of blades)
  • rotorsamplepoints_z::Array{TF,N}`: vertical wind location of points to sample across the rotor swept area when calculating the effective wind speed for the wind turbine. Points are centered at the hub (0,0) and scaled by the radius (1=tip of blades)
  • wind_resource::DiscretizedWindResource: wind resource discreption (directions, speeds, frequencies, etc)
  • model_set::AbstractModelSet: defines wake-realated models to be used in analysis
  • wind_farm_state_id::Int: index to correct state to use from wind resource provided. Defaults to 1
source
FLOWFarm.wake_count_iecMethod
wake_count_iec(turbinex, turbiney, winddirection, diameter; return_turbines=true)

Adapted from NREL's floris

Finds the number of turbines waking each turbine for the given
wind direction. Waked directions are determined using the formula
in Figure A.1 in Annex A of the IEC 61400-12-1:2017 standard.

Arguments

  • turbinex::Array{T,1}: x locations of turbines in global reference frame
  • turbiney::Array{T,1}: y locations of turbines in global reference frame
  • winddirection::Float: wind direction in radians in meteorological coordinates (0 rad. = from North)
  • diameter::Array{T,1}: diameters of all wind turbines
source
FLOWFarm.wake_deficit_modelMethod
wake_deficit_model(locx, locy, locz, turbine_x, turbine_y, turbine_z, deflection_y, deflection_z, upstream_turbine_id, downstream_turbine_id, hub_height, rotor_diameter, turbine_ai, turbine_local_ti, turbine_ct, turbine_yaw, model::GaussOriginal)

Computes the wake deficit at a given location using the Gaussian wake model presented by Bastankhah and Porte-Agel in the paper: "A new analytical model for wind-turbine wakes" (2014)

Arguments

  • locx::Float: x coordinate where wind speed is calculated
  • locy::Float: y coordinate where wind speed is calculated
  • locz::Float: z coordinate where wind speed is calculated
  • turbine_x::Array(Float): vector containing x coordinates for all turbines in farm
  • turbine_y::Array(Float): vector containing y coordinates for all turbines in farm
  • turbine_z::Array(Float): vector containing z coordinates for all turbines in farm
  • deflection_y::Float: deflection in the y direction of downstream wake
  • deflection_z::Float: deflection in the z direction of downstream wake
  • upstream_turbine_id::Int: index of the upstream wind turbine creating the wake
  • downstream_turbine_id::Int: index of the downstream turbine feeling the wake (if not referencing a turbine set to zero)
  • hub_height::Array(Float): vector containing hub heights for all turbines in farm
  • rotor_diameter::Array(Float): vector containing rotor diameters for all turbines in farm
  • turbine_ai::Array(Float): vector containing initial velocity deficits for all turbines in farm
  • turbine_local_ti::Array(Float): vector containing local turbulence intensities for all turbines in farm
  • turbine_ct::Array(Float): vector containing thrust coefficients for all turbines in farm
  • turbine_yaw::Array(Float): vector containing the yaw angle? for all turbines in farm
  • model::GaussOriginal: indicates the wake model in use
source
FLOWFarm.wake_deficit_modelMethod
wake_deficit_model(locx, locy, locz, turbine_x, turbine_y, turbine_z, deflection_y, deflection_z, upstream_turbine_id, downstream_turbine_id, hub_height, rotor_diameter, turbine_ai, turbine_local_ti, turbine_ct, turbine_yaw, model::GaussSimple)

Computes the wake deficit at a given location using the Gaussian wake model presented by Bastankhah and Porte-Agel in the paper: "A new analytical model for wind-turbine wakes" (2014) as modified for IEA Task 37 Case Studies 3 and 4

Arguments

  • locx::Float: x coordinate where wind speed is calculated
  • locy::Float: y coordinate where wind speed is calculated
  • locz::Float: z coordinate where wind speed is calculated
  • turbine_x::Array(Float): vector containing x coordinates for all turbines in farm
  • turbine_y::Array(Float): vector containing y coordinates for all turbines in farm
  • turbine_z::Array(Float): vector containing z coordinates for all turbines in farm
  • deflection_y::Float: deflection in the y direction of downstream wake
  • deflection_z::Float: deflection in the z direction of downstream wake
  • upstream_turbine_id::Int: index of the upstream wind turbine creating the wake
  • downstream_turbine_id::Int: index of the downstream turbine feeling the wake (if not referencing a turbine set to zero)
  • hub_height::Array(Float): vector containing hub heights for all turbines in farm
  • rotor_diameter::Array(Float): vector containing rotor diameters for all turbines in farm
  • turbine_ai::Array(Float): vector containing initial velocity deficits for all turbines in farm
  • turbine_local_ti::Array(Float): vector containing local turbulence intensities for all turbines in farm
  • turbine_ct::Array(Float): vector containing thrust coefficients for all turbines in farm
  • turbine_yaw::Array(Float): vector containing the yaw angle? for all turbines in farm
  • model::GaussSimple: indicates the wake model in use
source
FLOWFarm.wake_deficit_modelMethod
wake_deficit_model(locx, locy, locz, turbine_x, turbine_y, turbine_z, deflection_y, deflection_z, upstream_turbine_id, downstream_turbine_id, hub_height, rotor_diameter, turbine_ai, turbine_local_ti, turbine_ct, turbine_yaw, model::GaussYawVariableSpread)

Computes the wake deficit at a given location using the The Gaussian wake model presented by Bastankhah and Porte-Agel in the paper: "Experimental and theoretical study of wind turbine wakes in yawed conditions" (2016) The spread rate is adjusted based on local turbulence intensity as in Niayifar and Porte-Agel 2016

Arguments

  • locx::Float: x coordinate where wind speed is calculated
  • locy::Float: y coordinate where wind speed is calculated
  • locz::Float: z coordinate where wind speed is calculated
  • turbine_x::Array(Float): vector containing x coordinates for all turbines in farm
  • turbine_y::Array(Float): vector containing y coordinates for all turbines in farm
  • turbine_z::Array(Float): vector containing z coordinates for all turbines in farm
  • deflection_y::Float: deflection in the y direction of downstream wake
  • deflection_z::Float: deflection in the z direction of downstream wake
  • upstream_turbine_id::Int: index of the upstream wind turbine creating the wake
  • downstream_turbine_id::Int: index of the downstream turbine feeling the wake (if not referencing a turbine set to zero)
  • hub_height::Array(Float): vector containing hub heights for all turbines in farm
  • rotor_diameter::Array(Float): vector containing rotor diameters for all turbines in farm
  • turbine_ai::Array(Float): vector containing initial velocity deficits for all turbines in farm
  • turbine_local_ti::Array(Float): vector containing local turbulence intensities for all turbines in farm
  • turbine_ct::Array(Float): vector containing thrust coefficients for all turbines in farm
  • turbine_yaw::Array(Float): vector containing the yaw angle? for all turbines in farm
  • model::GaussYawVariableSpread: indicates the wake model in use
source
FLOWFarm.wake_deficit_modelMethod
wake_deficit_model(locx, locy, locz, turbine_x, turbine_y, turbine_z, deflection_y, deflection_z, upstream_turbine_id, downstream_turbine_id, hub_height, rotor_diameter, turbine_ai, turbine_local_ti, turbine_ct, turbine_yaw, model::GaussYaw)

Computes the wake deficit at a given location using the The Gaussian wake model presented by Bastankhah and Porte-Agel in the paper: "Experimental and theoretical study of wind turbine wakes in yawed conditions" (2016)

Arguments

  • locx::Float: x coordinate where wind speed is calculated
  • locy::Float: y coordinate where wind speed is calculated
  • locz::Float: z coordinate where wind speed is calculated
  • turbine_x::Array(Float): vector containing x coordinates for all turbines in farm
  • turbine_y::Array(Float): vector containing y coordinates for all turbines in farm
  • turbine_z::Array(Float): vector containing z coordinates for all turbines in farm
  • deflection_y::Float: deflection in the y direction of downstream wake
  • deflection_z::Float: deflection in the z direction of downstream wake
  • upstream_turbine_id::Int: index of the upstream wind turbine creating the wake
  • downstream_turbine_id::Int: index of the downstream turbine feeling the wake (if not referencing a turbine set to zero)
  • hub_height::Array(Float): vector containing hub heights for all turbines in farm
  • rotor_diameter::Array(Float): vector containing rotor diameters for all turbines in farm
  • turbine_ai::Array(Float): vector containing initial velocity deficits for all turbines in farm
  • turbine_local_ti::Array(Float): vector containing local turbulence intensities for all turbines in farm
  • turbine_ct::Array(Float): vector containing thrust coefficients for all turbines in farm
  • turbine_yaw::Array(Float): vector containing the yaw angle? for all turbines in farm
  • model::GaussYaw: indicates the wake model in use
source
FLOWFarm.wake_deficit_modelMethod
wake_deficit_model(locx, locy, locz, turbine_x, turbine_y, turbine_z, deflection_y, deflection_z, upstream_turbine_id, downstream_turbine_id, hub_height, rotor_diameter, turbine_ai, turbine_local_ti, turbine_ct, turbine_yaw, model::JensenCosine)

Computes the wake deficit according to the original Jensen cosine wake model, from the paper: "A Note on Wind Generator Interaction" by N.O. Jensen (1983)

Arguments

  • locx::Float: x coordinate where wind speed is calculated
  • locy::Float: y coordinate where wind speed is calculated
  • locz::Float: z coordinate where wind speed is calculated
  • turbine_x::Array(Float): vector containing x coordinates for all turbines in farm
  • turbine_y::Array(Float): vector containing y coordinates for all turbines in farm
  • turbine_z::Array(Float): vector containing z coordinates for all turbines in farm
  • deflection_y::Float: deflection in the y direction of downstream wake
  • deflection_z::Float: deflection in the z direction of downstream wake
  • upstream_turbine_id::Int: index of the upstream wind turbine creating the wake
  • downstream_turbine_id::Int: index of the downstream turbine feeling the wake (if not referencing a turbine set to zero)
  • hub_height::Array(Float): vector containing hub heights for all turbines in farm
  • rotor_diameter::Array(Float): vector containing rotor diameters for all turbines in farm
  • turbine_ai::Array(Float): vector containing initial velocity deficits for all turbines in farm
  • turbine_local_ti::Array(Float): vector containing local turbulence intensities for all turbines in farm
  • turbine_ct::Array(Float): vector containing thrust coefficients for all turbines in farm
  • turbine_yaw::Array(Float): vector containing the yaw angle? for all turbines in farm
  • model::JensenCosine: indicates the wake model in use
source
FLOWFarm.wake_deficit_modelMethod
wake_deficit_model(locx, locy, locz, turbine_x, turbine_y, turbine_z, deflection_y, deflection_z, upstream_turbine_id, downstream_turbine_id, hub_height, rotor_diameter, turbine_ai, turbine_local_ti, turbine_ct, turbine_yaw, model::JensenTopHat)

Computes the wake deficit according to the original Jensen top hat wake model, from the paper: "A Note on Wind Generator Interaction" by N.O. Jensen (1983)

Arguments

  • locx::Float: x coordinate where wind speed is calculated
  • locy::Float: y coordinate where wind speed is calculated
  • locz::Float: z coordinate where wind speed is calculated
  • turbine_x::Array(Float): vector containing x coordinates for all turbines in farm
  • turbine_y::Array(Float): vector containing y coordinates for all turbines in farm
  • turbine_z::Array(Float): vector containing z coordinates for all turbines in farm
  • deflection_y::Float: deflection in the y direction of downstream wake
  • deflection_z::Float: deflection in the z direction of downstream wake
  • upstream_turbine_id::Int: index of the upstream wind turbine creating the wake
  • downstream_turbine_id::Int: index of the downstream turbine feeling the wake (if not referencing a turbine set to zero)
  • hub_height::Array(Float): vector containing hub heights for all turbines in farm
  • rotor_diameter::Array(Float): vector containing rotor diameters for all turbines in farm
  • turbine_ai::Array(Float): vector containing initial velocity deficits for all turbines in farm
  • turbine_local_ti::Array(Float): vector containing local turbulence intensities for all turbines in farm
  • turbine_ct::Array(Float): vector containing thrust coefficients for all turbines in farm
  • turbine_yaw::Array(Float): vector containing the yaw angle? for all turbines in farm
  • model::JensenTopHat: indicates the wake model in use
source
FLOWFarm.wake_deficit_modelMethod
wake_deficit_model(locx, locy, locz, turbine_x, turbine_y, turbine_z, deflection_y, deflection_z, upstream_turbine_id, downstream_turbine_id, hub_height, rotor_diameter, turbine_ai, turbine_local_ti, turbine_ct, turbine_yaw, model::MultiZone)

Computes the wake deficit at a given location using the original MultiZone "FLORIS" wake model, from the paper: "Wind plant power optimization through yaw control using a parametric model for wake effects—a CFD simulation study" by Gebraad et al. (2014)

Arguments

  • locx::Float: x coordinate where wind speed is calculated
  • locy::Float: y coordinate where wind speed is calculated
  • locz::Float: z coordinate where wind speed is calculated
  • turbine_x::Array(Float): vector containing x coordinates for all turbines in farm
  • turbine_y::Array(Float): vector containing y coordinates for all turbines in farm
  • turbine_z::Array(Float): vector containing z coordinates for all turbines in farm
  • deflection_y::Float: deflection in the y direction of downstream wake
  • deflection_z::Float: deflection in the z direction of downstream wake
  • upstream_turbine_id::Int: index of the upstream wind turbine creating the wake
  • downstream_turbine_id::Int: index of the downstream turbine feeling the wake (if not referencing a turbine set to zero)
  • hub_height::Array(Float): vector containing hub heights for all turbines in farm
  • rotor_diameter::Array(Float): vector containing rotor diameters for all turbines in farm
  • turbine_ai::Array(Float): vector containing initial velocity deficits for all turbines in farm
  • turbine_local_ti::Array(Float): vector containing local turbulence intensities for all turbines in farm
  • turbine_ct::Array(Float): vector containing thrust coefficients for all turbines in farm
  • turbine_yaw::Array(Float): vector containing the yaw angle? for all turbines in farm
  • model::MultiZone: indicates the wake model in use
source
FLOWFarm.wake_deflection_modelMethod
wake_deflection_model(locx, locy, locz, turbine_x, turbine_yaw, turbine_ct, turbine_id, rotor_diameter, turbine_local_ti, model::GaussYawDeflection)

Calculates the horizontal deflection of the wind turbine wake

Based on:
[1] Bastankhah and Porte-Agel 2016 "Experimental and theoretical study of
wind turbine wakes in yawed conditions"
source
FLOWFarm.wake_deflection_modelMethod
wake_deflection_model(oc, turbine_x, turbine_yaw, turbine_ct, turbine_id, rotor_diameter, turbine_local_ti, model::GaussYawVariableSpreadDeflection)

Calculates the horizontal deflection of the wind turbine wake. Varies based on local turbulence intensity.

Based on:
[1] Bastankhah and Porte-Agel 2016 "Experimental and theoretical study of
wind turbine wakes in yawed conditions"
[2] Niayifar and Porte-Agel 2016 "Analytical Modeling of Wind Farms:
A New Approach for Power Prediction"
source
FLOWFarm.wake_deflection_modelMethod
wake_deflection_model(locx, locy, locz, turbine_id, turbine_definition::TurbineDefinition, model::JiminezYawDeflection)

Calculates the horizontal deflection of the wind turbine wake

Based on:
[1] Jiminez 2010 "Wake defl ection of a wind turbine in yaw"
[2] Gebraad 2014 "Wind plant optimization by yaw control using a parametric wake model"
this version ignores the corrections made to the yaw model for rotor rotation as described in [2] and
[3] Thomas 2017 "Improving the FLORIS wind plant model for compatibility with gradient-based optimization"
source
FLOWFarm.wake_deflection_modelMethod
wake_deflection_model(locx, locy, locz, turbine_id, turbine_definition::TurbineDefinition, model::NoYawDeflection, windfarmstate::SingleWindFarmState)

Bypasses yaw deflection calculations.
source
FLOWFarm.write_turb_loc_YAMLMethod
write_turb_loc_YAML(file_name, data)

write turbine locations and related information to .yaml

Arguments

  • file_name::String: path/and/name/of/location/file.yaml
source