Reference
FLOWFarm.CumulativeCurl
— TypeCumulativeCurl(a_s, b_s, c_s1, c_s2, a_f, b_f, c_f)
Container for parameters related to the Cumulative Curl model used in FLORIS v3 as shown in Bay 2022 (https://doi.org/10.5194/wes-2022-17)
Arguments
a_s::Float
: parameter relating turbulence intensity to wake spread. Default value is 0.179367259b_s::Float
: parameter controlling the default wake spread. Default value is 0.0118889215c_s1::Float
: parameter relating Ct to initial wake width. Default value is 0.0563691592c_s2::Float
: parameter affecting initial wake width. Default value is 0.13290157a_f::Float
: Default value is 3.11b_f::Float
: Default value is -0.68c_f::Float
: Default value is 2.41- 'wec_factor': paramter for the wake expansion continuation Default is [1.0]
FLOWFarm.DiscretizedWindResource
— TypeDiscritizedWindResource(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 radianswind_speeds::Array{Float,1}(Nstates)
: an array of wind speeds corresponding to each wind farm state in meters/secondwind_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 stateair_density::Float
: the air densityambient_ti::Array{Float,1}
: an array of the ambient turbulence intensity for each wind directionwind_shear_model::Array{AbstractWindShearModel}(1)
: contains a struct defining the desired turbulence intensity model
FLOWFarm.GaussOriginal
— TypeGaussOriginal(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
FLOWFarm.GaussSimple
— TypeGaussSimple(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 wakewec_factor::Array{Float}
: parameter artificial wake spreading for wake expansion continuation (WEC) optimization
FLOWFarm.GaussYaw
— TypeGaussYaw(turbulence_intensity, horizontal_spread_rate, vertical_spread_rate, alpha_star, beta_star, interpolation)
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.interpolation::Bool
: boolean stating if the the near wake should be interpolated. Default value is true.
FLOWFarm.GaussYawDeflection
— TypeGaussYawDeflection(horizontal_spread_rate, vertical_spread_rate, alpha_star, beta_star, interpolation)
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.interpolation::Bool
: boolean stating if the the near wake should be interpolated. Default value is true.
FLOWFarm.GaussYawVariableSpread
— TypeGaussYawVariableSpread(turbulence_intensity, horizontal_spread_rate, vertical_spread_rate, alpha_star, beta_star, interpolation)
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.interpolation::Bool
: boolean stating if the the near wake should be interpolated. Default value is true.
FLOWFarm.GaussYawVariableSpreadDeflection
— TypeGaussYawDeflectionVariableSpread(alpha_star, beta_star, k1, k2, interpolation)
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 intensityk2::Float
: second parameter tuning wake spread as based on turbulence intensityinterpolation::Bool
: boolean stating if the the near wake should be interpolated. Default value is true.
FLOWFarm.JensenCosine
— TypeJensenCosine(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.1beta::Float
: parameter controlling the width of the cosine function. Default value is 20.0 deg., given in radians.
FLOWFarm.JensenTopHat
— TypeJensenTopHat(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
FLOWFarm.JiminezYawDeflection
— TypeJiminezYawDeflection(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
FLOWFarm.Levelized
— TypeLevelized(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 moduleBOS::Float
: Balance of System (Costs outside of turbine i.e. operation and maintenance)FC::Float
: Financial Costs including construction and contingencyFCR::Float
: Fixed Charge RateOpEx::Float
: Operational Expenditures
FLOWFarm.LocalTIModelGaussTI
— TypeLocalTIModelGaussTI()
Calculate local turbulence intensity using the model presented in Qian and Ishihara (2018)
FLOWFarm.LocalTIModelMaxTI
— TypeLocalTIModelMaxTI(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 modelbstar::Float
: wake spreading parameter from Bastankhah and Porte-Agel Gaussian wake modelk1::Float
: slope of k vs TI curvek2::Float
: vertical offset of k vs TI curve
FLOWFarm.LocalTIModelNoLocalTI
— TypeLocalTIModelNoLocalTI()
Don't calculate local turbulence intensity. Ambient TI will be used instead for all points
FLOWFarm.MultiZone
— TypeMultiZone(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.065ke::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.
FLOWFarm.MultizoneDeflection
— TypeMultizoneDeflection(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.1ad::Float
:Helps define the horizontal deflection of the wake at 0 deg yawbd::Float
:Helps define the horizontal deflection of the wake due to downwind distance at 0 deg yaw
FLOWFarm.NoWakeDeficit
— TypeNoWakeDeficit()
FLOWFarm.NoYawDeflection
— TypeNoYawDeflection()
Allows for bypassing deflection calculations.
FLOWFarm.PowerLawWindShear
— TypePowerLawWindShear(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 shearground_height::Float
: height of the ground (typically zero)shear_order::Bool
: when shear should be calculated. Can be "first", "last", or "none"
FLOWFarm.PowerModelConstantCp
— TypePowerModelConstantCp(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
FLOWFarm.PowerModelCpPoints
— TypePowerModelCpPoints(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/scp_points::Array{N,Float}
: power coefficient values corresponding to the provided speeds- 'pp::TF': exponent for adjusting power for wind turbine yaw
FLOWFarm.PowerModelPowerCurveCubic
— TypePowerModelPowerCurveCubic()
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
FLOWFarm.PowerModelPowerPoints
— TypePowerModelPowerPoints(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/spower_points::Array{N,Float}
: power values corresponding to the provided speeds- 'pp::TF': exponent for adjusting power for wind turbine yaw
FLOWFarm.ThrustModelConstantCt
— TypeThrustModelConstantCt(ct::Float)
Stores a constant ct value for wake calculations
Arguments
ct::Float
: a constant ct value for computation
FLOWFarm.ThrustModelCtPoints
— TypeThrustModelCtPoints(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 turbinethrust_model::ThrustModelCtPoints
: Struct containing ct and velocity points for ct curve
FLOWFarm.WindFarm
— TypeWindFarm(windfarm, windresource, windfarmstates)
Struct defining a wind farm
Arguments
turbine_x::Array{Float}(Nturbines)
: contains windturbine x coordinates in the global reference frameturbine_y::Array{Float}(Nturbines)
: contains windturbine y coordinates in the global reference frameturbine_z::Array{Float}(Nturbines)
: contains windturbine base/z coordinates in the global reference frameturbine_definition_ids::Array{Int}(Nturbines)
: contains integers for each wind turbine specifying its definitionturbine_definitions::Array{AbstractTurbineDefinition}(Ntypes)
: contains structs defining each wind turbine definition (design) used in the farm
FLOWFarm.WindFarmModelSet
— TypeWindFarmModelSet(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 modelwake_deflection_model::AbstractWakeDeflectionModel
: contains a struct defining the desired wake deflection modelwake_combination_model::AbstractWakeCombinationModel
: contains a struct defining the desired wake combination modellocal_ti_model::AbstractTurbulenceIntensityModel
: contains a struct defining the desired turbulence intensity modelpoint_velocity_average_factor::Float
: factor used to determine how point velocity is averaged across the turbine (3 would result in a cubic mean)
FLOWFarm.boundary_struct
— Typeboundary_struct
Struct defining the boundary constraints
Arguments
turbine_x
: Vector containing x positions of turbinesturbine_y
: Vector containing y positions of turbinesboundary_scaling_factor
: Single float that scales the constraintboundary_function
: function that takes the boundary vector and the design variables to update the boundary vectorboundary_vec
: Vector containing the boundary constraintsjacobian
: Matrix containing the jacobian of the boundary constraintsconfig
: The ForwardDiff config object if using ForwardDiff for jacboian calculationupdate_function
: function that takes the design variables x and updates the boundary struct
FLOWFarm.preallocations_struct
— Typepreallocations_struct
struct that holds all the preallocated space for AEP calculation with one per thread used
Arguments
prealloc_turbine_velocities
: Vector containing preallocated space for turbine velocitiesprealloc_turbine_ct
: Vector containing preallocated space for turbine ctprealloc_turbine_ai
: Vector containing preallocated space for turbine aiprealloc_turbine_local_ti
: Vector containing preallocated space for turbine local tiprealloc_wake_deficits
: Matrix containing preallocated space for wake deficitsprealloc_contribution_matrix
: Matrix containing preallocated space for contribution matrixprealloc_deflections
: Matrix containing preallocated space for deflectionsprealloc_sigma_squared
: Matrix containing preallocated space for sigma squared
FLOWFarm.spacing_struct
— Typespacing_struct
Struct defining the spacing constraints
Arguments
turbine_x
: Vector containing x positions of turbinesturbine_y
: Vector containing y positions of turbinesconstraint_spacing
: Single float that defines the minimum spacing between turbines in metersconstraint_scaling
: Single float that scales the constraintspacing_vec
: Vector containing the spacing constraintsjacobian
: Matrix containing the jacobian of the spacing constraintsconfig
: The ForwardDiff config object if using ForwardDiff for jacboian calculationupdate_function
: function that takes the design variables x and updates the spacing struct
FLOWFarm.sparse_AEP_struct_stable_pattern
— TypesparseAEPstructstablepattern
Struct that holds all the necessary variables to calculate the AEP gradient using a stable sparsity pattern
Arguments
caches
: vector of SparseDiffTools jacobian caches for each wind statejacobians
: vector of sparse matracies containing jacobians for each wind statestate_gradients
: 2d array, each row is a state gradient (used for threads)turbine_powers
: 2d array that holds the powers or each turbine (used for threads)adtype
: AutoSparseForwardDiff object needed for SparseDiffTools
FLOWFarm.sparse_AEP_struct_unstable_pattern
— TypesparseAEPstructunstablepattern
Struct that holds all the necessary variables to calculate the AEP gradient using an unstable sparsity pattern
Arguments
deficit_thresholds
: Vector of floats that define the deficit thresholds for each wind statepatterns
: 3d array that holds the sparsity patterns for each wind statestate_gradients
: 2d array, each row is a state gradient (used for threads)jacobians
: Vector of sparse arrays containing jacobians for each wind stateturbine_powers
: 2d array that holds the powers or each turbine (used for threads)farm
: WindFarm structold_patterns
: 3d array that holds the old sparsity patterns for each wind statecolors
: 2d array that holds the colors for each wind statestate_powers
: 1d array that holds the state powerschunksize
: Chunksize for the AutoSparseForwardDiff object
FLOWFarm.sparse_boundary_struct
— Typesparseboundarystruct
Struct that holds all the necessary variables to calculate the boundary constraints using sparse methods
Arguments
turbine_x
: Vector containing x positions of turbinesturbine_y
: Vector containing y positions of turbinesjacobian
: Sparse matrix containing the jacobian of the boundary constraintsad
: AutoSparseForwardDiff objectcache
: SparseJacobianCache object for SparseDiffToolsboundary_vec
: Vector containing the boundary constraintsboundary_function
: Function that calculates the boundary constraintsupdate_function
: Function that updates the boundary struct with the new design variablesboundary_scaling_factor
: Single float that scales the boundary constraint
FLOWFarm.sparse_spacing_struct
— Typesparsespacingstruct
Struct that holds all the necessary variables to calculate the spacing constraints using sparse methods
Arguments
turbine_x
: Vector containing x positions of turbinesturbine_y
: Vector containing y positions of turbinesconstraint_spacing
: Single float that defines the minimum spacing between turbinesconstraint_scaling
: Single float that scales the constraintspacing_vec
: Vector containing the spacing constraintsjacobian
: Sparse matrix containing the jacobian of the spacing constraintscache
: SparseJacobianCache object for SparseDiffToolsupdate_function
: Function that updates the spacing struct with the new design variablesrelevant_list
: 2d array that holds the relevant turbine pairs for the spacing constraint (column 1 holds the first turbine and column 2 holds the second turbine in the pair)ad
: AutoSparseForwardDiff objectsafe_design_variables
: Vector containing the last set of design variables that satisfy the constraintsfull_spacing_vec
: Vector containing the full spacing constraints of the farm for final evaluation
FLOWFarm.wind_farm_constants_struct
— Typewindfarmconstants_struct
struct that holds all the constants for the wind farm
Arguments
turbine_z
: Vector containing z positions of ground the turbines sit onct_models
: Vector containing ct_models for each turbinegenerator_efficency
: Vector containing the generator efficiency of each turbinecut_in_speed
: Vector containing the cut in speed of each turbinecut_out_speed
: Vector containing the cut out speed of each turbinerated_speed
: Vector containing the rated speed of each turbinerated_power
: Vector containing the rated power of each turbinewind_resource
: The windresource structpower_models
: Vector containing power models of each turbinemodel_set
: The models_set structrotor_sample_points_y
: Vector containing y sample pointsrotor_sample_points_z
: Vector containing z sample points
FLOWFarm.wind_farm_struct
— Typewindfarmstruct
Unifying struct defining a wind farm and all necessary variables to calculate the AEP
Arguments
turbine_x
: Vector containing x positions of turbinesturbine_y
: Vector containing y positions of turbineshub_height
: Vector containing hub heights of each turbines as measured form the ground the turbines sit onturbine_yaw
: Vector containing yaw angle of each turbine in radiansrotor_diameter
: Vector containing the rotor diameter of each turbineresults
: DiffResults object to extract AEP when calculating AEP gradientconstants
: windfarmconstants_structAEP_scale
: Scaling factor for the AEPideal_AEP
: The ideal AEP of the farmpreallocations
: preallocated spaceupdate_function
: function that takes the design variables x and updates the farm structAEP_gradient
: The gradient of the AEPAEP
: The AEP of the farmconfig
: The ForwardDiff config object if using ForwardDiff for AEP gradient calculation, otherwise nothingforce_single_thread
: Boolean that forces the code to run in a single thread
FLOWFarm.DiscreteCircum
— FunctionPointsOnCircum(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 circlecenter_y::Float64
: cartesian y-coordinate for the center of the circler::Float64
: distance from circle's center to the circumference pointsn::Float64
: defaults to 100, is the number of discrete evenly-spaced points that will be returned along the circle's circumference
FLOWFarm.GaussianTI
— MethodGaussianTI(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. frameturbine_x::Array{Float,nTurbines}
: turbine wind direction locations in the wind direction reference frameturbine_y::Array{Float,nTurbines}
: turbine cross wind locations in the wind direction reference framerotor_diameter::Array{Float,nTurbines}
: rotor diameters of all turbineshub_height::Array{Float,nTurbines}
: hub heights of all turbines relative to the groundturbine_ct::Array{Float,nTurbines}
: thrust coefficient of each turbine for the given statesorted_turbine_index::Array{Float,nTurbines}
: turbine north-south locations in the global reference frameambient_ti::Float
: ambient turbulence intensitydiv_sigma::Float
: ?div_ti::Float
: ?
FLOWFarm.VR_boundary
— MethodVR_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 either 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 boundarybndry_y::Array{Float,1}
: 1-D array of y-coordinates for the vertices around a singlar closed boundarystart_dist::Float64
: the distance (positive or negative) along the boundary from the first boundary point where the turbines will begin to be placedturb_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.
FLOWFarm.VR_boundary_startup
— MethodVR_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 boundarybndry_y::Array{Float,1}
: 1-D array of y-coordinates for the vertices around a singlar closed boundarystart_dist::Float64
: the distance (positive or negative) along the boundary from the first boundary point where the turbines will begin to be placedturb_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
FLOWFarm._ct_to_axial_ind_func
— Method_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
FLOWFarm._gauss_yaw_potential_core
— Method_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.
FLOWFarm._gauss_yaw_spread
— Method_gauss_yaw_spread(dt, k, dx, x0, yaw)
Helper function for wakedeficitmodel when using the GaussYaw model. Computes the standard deviation of the wake.
FLOWFarm._gauss_yaw_spread_interpolated
— Method_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.
FLOWFarm._niayifar_added_ti_function
— Method_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 interestd_dst::Float
: downstream turbine rotor diameterd_ust::Float
: upstream turbine rotor diameterh_ust::Float
: upstream turbine hub heighth_dst::Float
: downstream turbine hub heightct_ust::Float
: upstream turbine thrust coefficientkstar_ust::Float
: upstream turbine wake expansion ratedelta_y::Float
: cross wind separation from turbine to point of interestti_amb::Float
: ambient turbulence intensityti_ust::Float
: upstream turbine local turbulence intensityti_dst::Float
: downstream turbine local turbulence intensityti_area_ratio_in::Float
: current value of TI-area ratio for use in calculatin local TIs::Float
: smooth max smootheness parameter
FLOWFarm._remove_out_of_bounds_points
— Method_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 locationsz::AbstractArray
: vertical point locationsuse_perimeter_points::Bool
: flag that determines whether or not to include points on the boundary of the rotor-swept area
FLOWFarm.adjust_for_wind_shear
— Methodadjust_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 velocityreference_velocity::Float
: known velocity at reference_heightreference_height::Float
: height of known velocityground_height::Float
: height of the ground (typically zero)model::AbstractWindShearModel
: wind shear model to use for calculations
FLOWFarm.boundary_normals_calculator
— Methodboundary_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, counterclockwisenboundaries::Int
: the number of boundaries in the set
FLOWFarm.build_boundary_struct
— Methodbuildboundarystruct
buildboundarystruct(x,nturbines,nconstraints,scaling,constraintfunction,updatefunction;using_sparsity=true)
Arguments
x
: Vector containing the desired design variables for optimizationn_turbines
: Number of turbines in the wind farmn_constraints
: Number of boundary constraints (n_turbines * number of sides of the boundary)scaling
: The scaling factor for the boundary constraintconstraint_function
: Function that calculates the boundary constraintsupdate_function
: Function that updates the boundary struct with the new design variablesusing_sparsity
: Boolean to use sparsity in the jacobian calculation (default is true)
FLOWFarm.build_relevant_list
— Methodbuildrelevantlist(turbinex,turbiney,space,factor)
Helper function that builds the relevant list for the sparse spacing constraints
Arguments
turbine_x
: Vector containing x positions of turbinesturbine_y
: Vector containing y positions of turbinesspace
: Single float that defines the minimum spacing between turbinesfactor
: Single float that defines the factor of space to be considered relevant
FLOWFarm.build_spacing_struct
— Methodbuildspacingstruct(x,nturbines,space,scale,updatefunction)
function to build a spacing_struct
Arguments
x
: Vector containing the desired design variables for optimizationn_turbines
: Number of turbines in the wind farmspace
: The minimum spacing between turbinesscale
: The scaling factor for the spacing constraintupdate_function
: Function that updates the spacing struct with the new design variables
FLOWFarm.build_sparse_spacing_struct
— Methodbuildsparsespacingstruct(x,turbinex,turbiney,space,scale,updatefunction;firstopt=true,relevantspacing_factor=2)
Function that builds a sparsespacingstruct
Arguments
x
: Vector containing the design variablesturbine_x
: Vector containing x positions of turbinesturbine_y
: Vector containing y positions of turbinesspace
: Single float that defines the minimum spacing between turbinesscale
: Single float that scales the constraintupdate_function
: Function that updates the spacing struct with the new design variablesfirst_opt
: Boolean to determine if this is the first optimization (if true uses no spacing constraints)relevant_spacing_factor
: Single float that defines the factor of space to be considered relevant
FLOWFarm.build_stable_sparse_struct
— Methodbuildstablesparsestruct(x,turbinex,turbiney,turbinez,hubheight,turbineyaw,rotordiameter, ctmodels,generatorefficiency,cutinspeed,cutoutspeed,ratedspeed,ratedpower,windresource, powermodels,modelset,updatefunction;rotorsamplepointsy=[0.0],rotorsamplepointsz=[0.0], AEPscale=0.0,optx=false,opty=false,opthub=false,optyaw=false,optdiam=false,tolerance=1E-16, forcesingle_thread=false)
Function that builds a windfarmstruct and a sparseAEPstructstablepattern struct that goes with it
Arguments
x
: Vector containing the design variablesturbine_x
: Vector containing x positions of turbinesturbine_y
: Vector containing y positions of turbinesturbine_z
: Vector containing z positions of turbines (this is where the base meets the ground)hub_height
: Vector containing hub heights of turbinesturbine_yaw
: Vector containing yaw angles of turbinesrotor_diameter
: Vector containing rotor diameters of turbinesct_models
: Vector containing ct_models for each turbinegenerator_efficiency
: Vector containing generator efficiencies for each turbinecut_in_speed
: Vector containing cut in speeds for each turbinecut_out_speed
: Vector containing cut out speeds for each turbinerated_speed
: Vector containing rated speeds for each turbinerated_power
: Vector containing rated powers for each turbinewind_resource
: The DiscretizedWindResource structpower_models
: Vector containing power models for each turbinemodel_set
: The WindFarmModelSet for the wind farmupdate_function
: Function that updates the wind farm struct with the new design variablesrotor_sample_points_y
: Vector containing the y positions of the rotor sample pointsrotor_sample_points_z
: Vector containing the z positions of the rotor sample pointsAEP_scale
: Single float that scales the AEP, if 0.0 will be set to 1.0/ideal_AEPinput_type
: default is nothing and will be set to the type of x, if "ForwardDiff" then the input type will be set to ForwardDiff.dualopt_x
: Boolean to optimize x positions of turbinesopt_y
: Boolean to optimize y positions of turbinesopt_hub
: Boolean to optimize hub heights of turbinesopt_yaw
: Boolean to optimize yaw angles of turbinesopt_diam
: Boolean to optimize rotor diameters of turbinestolerance
: Single float that defines the tolerance for the jacobian pattern (default is 1E-16), set to 0.0 to use traditional sparsityforce_single_thread
: Boolean to force single thread calculation
FLOWFarm.build_stable_sparse_struct
— Methodbuildstablesparse_struct(x,farm;tolerance=1E-16)
Helper function that builds a sparseAEPstructstablepattern struct
Arguments
x
: Vector containing the design variablesfarm
: WindFarm structtolerance
: Single float that defines the tolerance for the jacobian pattern
FLOWFarm.build_unstable_sparse_struct
— Methodbuildunstablesparsestruct(x,turbinex,turbiney,turbinez,hubheight,turbineyaw,rotordiameter, ctmodels,generatorefficiency,cutinspeed,cutoutspeed,ratedspeed,ratedpower,windresource, powermodels,modelset,updatefunction;rotorsamplepointsy=[0.0],rotorsamplepointsz=[0.0], AEPscale=0.0,optx=false,opty=false,opthub=false,optyaw=false,optdiam=false,tolerance=1E-16, forcesingle_thread=false)
Function that builds a windfarmstruct and a sparseAEPstructunstablepattern struct
Arguments
x
: Vector containing the design variablesturbine_x
: Vector containing x positions of turbinesturbine_y
: Vector containing y positions of turbinesturbine_z
: Vector containing z positions of turbines (this is where the base meets the ground)hub_height
: Vector containing hub heights of turbinesturbine_yaw
: Vector containing yaw angles of turbinesrotor_diameter
: Vector containing rotor diameters of turbinesct_models
: Vector containing ct_models for each turbinegenerator_efficiency
: Vector containing generator efficiencies for each turbinecut_in_speed
: Vector containing cut in speeds for each turbinecut_out_speed
: Vector containing cut out speeds for each turbinerated_speed
: Vector containing rated speeds for each turbinerated_power
: Vector containing rated powers for each turbinewind_resource
: The DiscretizedWindResource structpower_models
: Vector containing power models for each turbinemodel_set
: The WindFarmModelSet for the wind farmupdate_function
: Function that updates the wind farm struct with the new design variablesrotor_sample_points_y
: Vector containing the y positions of the rotor sample pointsrotor_sample_points_z
: Vector containing the z positions of the rotor sample pointsAEP_scale
: Single float that scales the AEP, if 0.0 will be set to 1.0/ideal_AEPinput_type
: default is nothing and will be set to the type of x, if "ForwardDiff" then the input type will be set to ForwardDiff.dualopt_x
: Boolean to optimize x positions of turbinesopt_y
: Boolean to optimize y positions of turbinesopt_hub
: Boolean to optimize hub heights of turbinesopt_yaw
: Boolean to optimize yaw angles of turbinesopt_diam
: Boolean to optimize rotor diameters of turbinestolerance
: Single float that defines the tolerance for the jacobian pattern (default is 1E-16), set to 0.0 to use traditional sparsityforce_single_thread
: Boolean to force single thread calculation
FLOWFarm.build_unstable_sparse_struct
— Methodbuildunstablesparsestruct(x,farm,farmforwarddiff;tolerance=1E-16)
Helper function that builds a sparseAEPstructunstablepattern struct
Arguments
x
: Vector containing the design variablesfarm
: WindFarm structfarm_forwarddiff
: WindFarm struct with ForwardDiff input type for deficit tolerance calculationtolerance
: Single float that defines the tolerance for the jacobian pattern
FLOWFarm.build_wind_farm_struct
— Methodbuildwindfarmstruct(x,turbinex,turbiney,turbinez,hubheight,turbineyaw,rotordiameter, ctmodels,generatorefficiency,cutinspeed,cutoutspeed,ratedspeed,ratedpower,windresource, powermodels,modelset,updatefunction;rotorsamplepointsy=[0.0],rotorsamplepointsz=[0.0], AEPscale=0.0,inputtype=nothing,optx=false,opty=false,opthub=false,optyaw=false,optdiam=false, forcesinglethread=false)
function to build a windfarmstruct
Arguments
x
: Vector containing the desired design variables for optimization (if any)turbine_x
: Vector containing x positions of turbinesturbine_y
: Vector containing y positions of turbinesturbine_z
: Vector containing z positions of turbines (this is where the base meets the ground)hub_height
: Vector containing hub heights of turbinesturbine_yaw
: Vector containing yaw angles of turbinesrotor_diameter
: Vector containing rotor diameters of turbinesct_models
: Vector containing ct_models for each turbinegenerator_efficiency
: Vector containing generator efficiencies for each turbinecut_in_speed
: Vector containing cut in speeds for each turbinecut_out_speed
: Vector containing cut out speeds for each turbinerated_speed
: Vector containing rated speeds for each turbinerated_power
: Vector containing rated powers for each turbinewind_resource
: The DiscretizedWindResource structpower_models
: Vector containing power models for each turbinemodel_set
: The WindFarmModelSet for the wind farmupdate_function
: Function that updates the wind farm struct with the new design variablesrotor_sample_points_y
: Vector containing the y positions of the rotor sample pointsrotor_sample_points_z
: Vector containing the z positions of the rotor sample pointsAEP_scale
: Single float that scales the AEP, if 0.0 will be set to 1.0/ideal_AEPinput_type
: default is nothing and will be set to the type of x, if "ForwardDiff" then the input type will be set to ForwardDiff.dualopt_x
: Boolean to optimize x positions of turbinesopt_y
: Boolean to optimize y positions of turbinesopt_hub
: Boolean to optimize hub heights of turbinesopt_yaw
: Boolean to optimize yaw angles of turbinesopt_diam
: Boolean to optimize rotor diameters of turbinesforce_single_thread
: Boolean to force single thread calculation
FLOWFarm.calcMinorAngle
— FunctioncalcMinorAngle(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 boundarybndry_y::Array{Float,1}
: 1-D array of y-coordinates for the vertices around a singlar closed boundarybndry_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
FLOWFarm.calcSmallestAngle
— MethodcalcSmallestAngle(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 boundarybndry_y::Array{Float,1}
: 1-D array of y-coordinates for the vertices around a singlar closed boundary
FLOWFarm.calc_moment_stress
— Functioncalc_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 momentmy::Float
: y momentdx::Float
: x distance to the location of interestdy::Float
: y distance to the location of interest
Keyword Arguments
Rcyl::Float
: radius of the cylindertcyl::Float
: thickenss of the cylinder
FLOWFarm.calculate_aep!
— Methodcalculate_aep!
function calculate_aep!(farm,x)
Arguments
farm
: The windfarmstructx
: Vector containing the design variables
FLOWFarm.calculate_aep
— Methodcalculate_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 frameturbine_y::Array{Float,nTurbines}
: turbine north-south locations in the global reference frameturbine_z::Array{Float,nTurbines}
: turbine base height in the global reference framerotor_diameter::Array{Float,nTurbines}
hub_height::Array{TF,nTurbines}
: turbine hub heightsturbine_yaw::Array{TF,nTurbines}
: turbine yaw for the given wind direction in radiansct_model::AbstractThrustCoefficientModel
: defines how the thrust coefficient changes with state etcgenerator_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 AbstractPowerModelmodel_set::AbstractModelSet
: defines wake-realated models to be used in analysisrotor_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)
FLOWFarm.calculate_aep_gradient!
— Methodcalculateaepgradient!(farm,x)
function to calculate the AEP and its gradient for the wind farm, results are stored within the windfarmstruct and take into account the scaling factor
Arguments
farm
: The windfarmstructx
: Vector containing the design variables
FLOWFarm.calculate_aep_gradient!
— Methodcalculateaepgradient!(farm,x,sparse_struct::T)
function to calculate the AEP and its gradient for the wind farm using sparse methods
Arguments
farm
: The windfarmstructx
: Vector containing the design variables
FLOWFarm.calculate_aep_gradient!
— Methodcalculateaepgradient!(farm,x,sparse_struct::T)
Function that calculates the AEP gradient using a stable sparsity pattern
Arguments
farm
: WindFarm structx
: Vector containing the scaled design variablessparse_struct
: sparseAEPstructstablepattern struct
FLOWFarm.calculate_aep_gradient!
— Methodcalculateaepgradient!(farm,x,sparse_struct::T)
Function that calculates the AEP gradient using an unstable sparsity pattern
Arguments
farm
: WindFarm structx
: Vector containing the scaled design variablessparse_struct
: sparseAEPstructunstablepattern struct
FLOWFarm.calculate_boundary!
— Methodcalculateboundary!(boundaryvec,x,boundary_struct)
function to calculate the boundary constraints for the wind farm
Arguments
boundary_vec
: Vector containing the boundary constraints (just us boundarystruct.boundaryvec)x
: Vector containing the design variablesboundary_struct
: The boundary_struct
FLOWFarm.calculate_boundary_jacobian!
— Methodcalculateboundaryjacobian!(boundary_struct,x)
function to calculate the boundary constraints and the jacobian for the wind farm, results stored in boundary_struct
Arguments
boundary_struct
: The boundary_structx
: Vector containing the design variables
FLOWFarm.calculate_boundary_jacobian!
— Methodcalculateboundaryjacobian!(boundary_struct::T,x)
Function that builds a sparseboundarystruct
Arguments
boundary_struct
: sparseboundarystructx
: Vector containing the design variables
FLOWFarm.calculate_ct
— Methodcalculate_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
FLOWFarm.calculate_ct
— Methodcalculate_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 turbinethrust_model::ThrustModelCtPoints
: Struct containing ct and velocity points for ct curve
FLOWFarm.calculate_flow_field
— Methodcalculateflowfield(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 frameyrange::Range
: range defining north-west locations to sample in global reference framezrange::Range
: range defining vertical locations to sample in global reference framemodel_set::AbstractModelSet
: defines wake-realated models to be used in analysisturbine_x::Array{TF,nTurbines}
: turbine east-west locations in the global reference frameturbine_y::Array{TF,nTurbines}
: turbine north-south locations in the global reference frameturbine_z::Array{TF,nTurbines}
: turbine base height in the global reference frameturbine_yaw::Array{TF,nTurbines}
: turbine yaw for the given wind direction in radiansturbine_ct::Array{TF,nTurbines}
: thrust coefficient of each turbine for the given stateturbine_ai::Array{TF,nTurbines}
: turbine axial induction for the given staterotor_diameter::Array{TF,nTurbines}
: turbine rotor diametershub_height::Array{TF,nTurbines}
: turbine hub heightsturbine_local_ti::Array{TF,nTurbines}
: turbine local turbulence intensity for the given statesorted_turbine_index::Array{TF,nTurbines}
: turbine north-south locations in the global reference framewtvelocities::Array{TF,nTurbines}
: effective inflow wind speed for given statewind_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
FLOWFarm.calculate_ideal_aep
— Methodcalculate_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 ideal wind farm AEP (AEP with no wake loss)
Arguments
turbine_x::Array{Float,nTurbines}
: turbine east-west locations in the global reference frameturbine_y::Array{Float,nTurbines}
: turbine north-south locations in the global reference frameturbine_z::Array{Float,nTurbines}
: turbine base height in the global reference framerotor_diameter::Array{Float,nTurbines}
hub_height::Array{TF,nTurbines}
: turbine hub heightsturbine_yaw::Array{TF,nTurbines}
: turbine yaw for the given wind direction in radiansct_model::AbstractThrustCoefficientModel
: defines how the thrust coefficient changes with state etcgenerator_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 AbstractPowerModelmodel_set::AbstractModelSet
: defines wake-realated models to be used in analysisrotor_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)
FLOWFarm.calculate_local_ti
— Methodcalculate_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::LocalTIModelGaussTI; turbine_id=1, tol=1E-6)
Returns local turbulence intensity calculated using methods in Qian 2018 from the Journal of Wind Energy https://doi.org/10.1016/j.jweia.2018.04.010 with modification to account for yaw coming from Qian 2018 from Energies doi:10.3390/en11030665
Arguments
turbine_x::Array{Float,nTurbines}
: turbine wind direction locations in the wind direction reference frameturbine_y::Array{Float,nTurbines}
: turbine cross wind locations in the wind direction reference frameambient_ti::Float
: ambient turbulence intensityrotor_diameter::Array{Float,nTurbines}
: rotor diameters of all turbineshub_height::Array{Float,nTurbines}
: hub heights of all turbines relative to the groundturbine_yaw::Array{Float,nTurbines}
: yaw of all turbines for the current wind state in radiansturbine_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 frameturbine_inflow_velcities::Array{Float,nTurbines}
: effective inflow wind speed at each turbine for given stateturbine_ct::Array{Float,nTurbines}
: thrust coefficient of each turbine for the given stateti_model::LocalTIModelGaussTI
: contains a struct defining the desired turbulence intensity modelturbine_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
FLOWFarm.calculate_local_ti
— Methodcalculate_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 frameturbine_y::Array{Float,nTurbines}
: turbine cross wind locations in the wind direction reference frameambient_ti::Float
: ambient turbulence intensityrotor_diameter::Array{Float,nTurbines}
: rotor diameters of all turbineshub_height::Array{Float,nTurbines}
: hub heights of all turbines relative to the groundturbine_yaw::Array{Float,nTurbines}
: yaw of all turbines for the current wind state in radiansturbine_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 frameturbine_inflow_velcities::Array{Float,nTurbines}
: effective inflow wind speed at each turbine for given stateturbine_ct::Array{Float,nTurbines}
: thrust coefficient of each turbine for the given stateti_model::LocalTIModelMaxTI
: contains a struct defining the desired turbulence intensity model, no local TI in this caseturbine_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
FLOWFarm.calculate_local_ti
— Methodcalculate_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 frameturbine_y::Array{Float,nTurbines}
: turbine cross wind locations in the wind direction reference frameambient_ti::Float
: ambient turbulence intensityrotor_diameter::Array{Float,nTurbines}
: rotor diameters of all turbineshub_height::Array{Float,nTurbines}
: hub heights of all turbines relative to the groundturbine_yaw::Array{Float,nTurbines}
: yaw of all turbines for the current wind state in radiansturbine_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 frameturbine_inflow_velcities::Array{Float,nTurbines}
: effective inflow wind speed at each turbine for given stateturbine_ct::Array{Float,nTurbines}
: thrust coefficient of each turbine for the given stateti_model::LocalTIModelNoLocalTI
: contains a struct defining the desired turbulence intensity model, no local TI in this caseturbine_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
FLOWFarm.calculate_power
— Methodcalculate_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 generatorair_density::Float
: Air densityrotor_area::Float
: Rotor-swept area of the wind turbinewt_velocity::Float
: Inflow velocity to the wind turbinecut_in_speed::Float
: cut in speed of the wind turbinerated_speed::Float
: rated speed of the wind turbinecut_out_speed::Float
: cut out speed of the wind turbinerated_power::Float
: rated power of the wind turbinepower_model::PowerModelConstantCp
: Struct containing the cp value to be used in region 2
FLOWFarm.calculate_power
— Methodcalculate_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 generatorair_density::Float
: Air densityrotor_area::Float
: Rotor-swept area of the wind turbinewt_velocity::Float
: Inflow velocity to the wind turbinecut_in_speed::Float
: cut in speed of the wind turbinerated_speed::Float
: rated speed of the wind turbinecut_out_speed::Float
: cut out speed of the wind turbinerated_power::Float
: rated power of the wind turbinepower_model::PowerModelCpPoints
: Struct containing the velocity and cp values defining the cp curve
FLOWFarm.calculate_power
— Methodcalculate_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 generatorair_density::Float
: Air densityrotor_area::Float
: Rotor-swept area of the wind turbinewt_velocity::Float
: Inflow velocity to the wind turbinecut_in_speed::Float
: cut in speed of the wind turbinerated_speed::Float
: rated speed of the wind turbinecut_out_speed::Float
: cut out speed of the wind turbinerated_power::Float
: rated power of the wind turbinepower_model::PowerModelPowerCurveCubic
: Empty struct
FLOWFarm.calculate_power
— Methodcalculate_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 generatorair_density::Float
: Air densityrotor_area::Float
: Rotor-swept area of the wind turbinewt_velocity::Float
: Inflow velocity to the wind turbinecut_in_speed::Float
: cut in speed of the wind turbinerated_speed::Float
: rated speed of the wind turbinecut_out_speed::Float
: cut out speed of the wind turbinerated_power::Float
: rated power of the wind turbinepower_model::PowerModelPowerPoints
: Struct containing the velocity and power values defining the power curve
FLOWFarm.calculate_power_from_cp
— Methodcalculate_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 generatorair_density::Float
: Air densityrotor_area::Float
: Rotor-swept area of the wind turbinecp::Float
: Power coefficient of the wind turbinewt_velocity::Float
: Inflow velocity to the wind turbine
FLOWFarm.calculate_spacing!
— Methodcalculatespacing!(spacingvec,x,spacing_struct)
function to calculate the spacing constraints for the wind farm
Arguments
spacing_vec
: Vector containing the spacing constraints (just us spacingstruct.spacingvec)x
: Vector containing the design variablesspacing_struct
: The spacing_struct
FLOWFarm.calculate_spacing!
— Methodbuildspacingstruct(x,nturbines,space,scale,updatefunction)
Calculates the spacing constraints using sparse methods
FLOWFarm.calculate_spacing_jacobian!
— Methodcalculatespacingjacobian!(spacing_struct,x)
function to calculate the spacing constraints and the jacobian for the wind farm, results stored in spacing_struct
Arguments
spacing_struct
: The spacing_structx
: Vector containing the design variables
FLOWFarm.calculate_spacing_jacobian!
— Methodcalculatespacingjacobian!(spacing_struct,x)
Function that calculates the spacing constraints jacobian using sparse methods
Arguments
spacing_struct
: sparsespacingstructx
: Vector containing the design variables
FLOWFarm.calculate_state_aeps
— Methodcalculate_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 frameturbine_y::Array{Float,nTurbines}
: turbine north-south locations in the global reference frameturbine_z::Array{Float,nTurbines}
: turbine base height in the global reference framerotor_diameter::Array{Float,nTurbines}
hub_height::Array{TF,nTurbines}
: turbine hub heightsturbine_yaw::Array{TF,nTurbines}
: turbine yaw for the given wind direction in radiansct_model::AbstractThrustCoefficientModel
: defines how the thrust coefficient changes with state etcgenerator_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 AbstractPowerModelmodel_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)
FLOWFarm.calculate_state_turbine_powers
— Methodcalculate_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 frameturbine_y::Array{Float,nTurbines}
: turbine north-south locations in the global reference frameturbine_z::Array{Float,nTurbines}
: turbine base height in the global reference framerotor_diameter::Array{Float,nTurbines}
hub_height::Array{TF,nTurbines}
: turbine hub heightsturbine_yaw::Array{TF,nTurbines}
: turbine yaw for the given wind direction in radiansct_model::AbstractThrustCoefficientModel
: defines how the thrust coefficient changes with state etcgenerator_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 AbstractPowerModelmodel_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)
FLOWFarm.calculate_threshold
— Methodcalculatethreshold(x,farmforwarddiff,farm,tolerance,pow,stateid;preallocid=1,lock=nothing)
Helper function that calculates the deficit threshold for a single wind state
Arguments
x
: Vector containing the design variablesfarm_forwarddiff
: WindFarm struct with ForwardDiff input type for deficit tolerance calculationfarm
: WindFarm structtolerance
: Single float that defines the tolerance for the jacobian patternpow
: 1d array that holds the powers or each turbinestate_id
: Wind state idprealloc_id
: Preallocation id (to select the correct preallocated memory inside the wind farm struct)lock
: SpinLock object to lock the farm struct for multithreadeding
FLOWFarm.calculate_thresholds!
— Methodcalculatethresholds!(jacobians,thresholds,x,farmforwarddiff,farm,tolerance,pow,n_states)
Helper function that calculates the thresholds for each wind state
Arguments
jacobians
: Vector of sparse arrays holding jacobiansthresholds
: Vector of floats that define the deficit thresholds for each wind statex
: Vector containing the design variablesfarm_forwarddiff
: WindFarm struct with ForwardDiff input type for deficit tolerance calculationfarm
: WindFarm structtolerance
: Single float that defines the tolerance for the jacobian patternpow
: 2d array that holds the powers or each turbine (used for threads)n_states
: Number of wind states
FLOWFarm.calculate_turbine_power
— Methodcalculate_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
FLOWFarm.calculate_unstable_sparse_jacobian!
— Methodcalculateunstablesparsejacobian!(sparsestruct::T,x,farm,windstateid,prealloc_id,lock)
Helper function that calculates the sparse jacobian for a single wind state using an unstable sparsity pattern
Arguments
sparse_struct
: sparseAEPstructunstablepattern structx
: Vector containing the design variablesfarm
: WindFarm structwind_state_id
: Wind state idprealloc_id
: Preallocation id (to select the correct preallocated memory inside the wind farm struct)lock
: SpinLock object to lock the farm struct for multithreadeding
FLOWFarm.calculate_unstable_sparsity_pattern!
— Methodcalculateunstablesparsitypattern!(sparsestruct::T,x,windstateid,prealloc_id)
Helper function that calculates the sparsity pattern for a single wind state using an unstable sparsity pattern
Arguments
sparse_struct
: sparseAEPstructunstablepattern structx
: Vector containing the design variableswind_state_id
: Wind state idprealloc_id
: Preallocation id (to select the correct preallocated memory inside the wind farm struct)
FLOWFarm.calculate_wind_state_power!
— Methodcalculatewindstatepower!(pow,x,farm,stateid;preallocid=1,hoursper_year=365.25*24.0,lock=nothing)
Helper function that calculates the power fora a single wind state
Arguments
pow
: 1d array that holds the powers or each turbinex
: Vector containing the design variablesfarm
: WindFarm structstate_id
: Wind state idprealloc_id
: Preallocation id (to select the correct preallocated memory inside the wind farm struct)hours_per_year
: Single float that defines the hours per yearlock
: SpinLock object to lock the farm struct for multithreadeding
FLOWFarm.circle_boundary!
— Methodcircle_boundary!(center,radius,turbine_x,turbine_y)
calculate the distance squared from each turbine to a circular boundary. Negative means the turbine is inside the boundary
Arguments
boundary_vec
: vector containing distances from the boundarycenter::Float
: circular boundary center [x,y]radius::Float
: circulat boundary radiusturbine_x::Array{Float}
: turbine x locationsturbine_y::Array{Float}
: turbine y locations
FLOWFarm.circle_boundary
— Methodcircle_boundary(center,radius,turbine_x,turbine_y)
calculate the distance squared 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 radiusturbine_x::Array{Float}
: turbine x locationsturbine_y::Array{Float}
: turbine y locations
FLOWFarm.closeBndryList
— MethodcloseBndryList(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 boundarybndryPts_y::Array{Float,1}
: 1-D array of y-coordinates for the vertices around a singlar closed boundary
FLOWFarm.closeBndryLists
— MethodcloseBndryLists(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 boundariesbndryPts_y::Array{Float,1}
: N-D array of y-coordinates for the vertices around N-many closed boundaries
FLOWFarm.coordDist
— MethodcoordDist(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 pointy1::Float64
: y-coord of the first pointx2::Float64
: x-coord of the second pointy2::Float64
: y-coord of the second point
FLOWFarm.define_patterns!
— Methoddefinepatterns!(jacobians,x,farm,tolerance,pow,nstates)
Helper function that defines the jacobian patterns for each wind state
Arguments
jacobians
: Vector of sparse arrays holding jacobiansx
: Vector containing the design variablesfarm
: WindFarm structtolerance
: Single float that defines the tolerance for the jacobian patternpow
: 2d array that holds the powers or each turbine (used for threads)n_states
: Number of wind states
FLOWFarm.define_stable_jacobian_pattern
— Methoddefinestablejacobianpattern(x,farm,tolerance,pow,stateid;prealloc_id=1)
Helper function that defines the jacobian pattern for a single wind state
Arguments
x
: Vector containing the design variablesfarm
: WindFarm structtolerance
: Single float that defines the tolerance for the jacobian patternpow
: 1d array that holds the powers or each turbinestate_id
: Wind state idprealloc_id
: Preallocation id (to select the correct preallocated memory inside the wind farm struct)
FLOWFarm.distributed_velocity_op
— Functiondistributed_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 rOmega::Float
: rotor rotational speedr::Array{Float}
: radial locations of interestprecone::Float
: rotor precone angleyaw::Float
: rotor yaw angletilt::Float
: rotor tilt angleazimuth::Float
: blade azimuth anglerho::Float
: air density
Keyword Arguments
mu::Float
: air viscocity (can usually use the default)asound::Float
: speed of sound (can usually use the default)
FLOWFarm.find_upstream_turbines
— Methodfind_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 frameturbiney::Array{T,1}
: y locations of turbines in global reference framewinddirection::Real
orwinddirection::AbstractArray
: wind direction in radians in meteorological coordinates (0 rad. = from North)diameter::Array{T,1}
: diameters of all wind turbines
FLOWFarm.find_xyz_simple
— Methodfind_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 huby_hub::Float
: y location of hubz_hub::Float
: z location of hub (hub height if no topology)r::Array{Float}
: radial locations of interestprecone::Float
: rotor precone angleyaw::Float
: rotor yaw angleazimuth::Float
: blade azimuth angle
FLOWFarm.getNextFileName
— FunctiongetNextFileName(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 desiredfile_type::String
: ex "yaml", "txt", "csv", etc...max_check::Int
: the maximum number of files to check
FLOWFarm.getPerimeterLength
— MethodgetPerimeterLength(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 boundarybndry_y::Array{Float,1}
: 1-D array of y-coordinates for the vertices around a singlar closed boundary
FLOWFarm.getUpDwnYvals
— MethodgetUpDwnYvals(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 examinedbndry_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"
FLOWFarm.get_boundary_yaml
— Methodget_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
FLOWFarm.get_moments
— Methodget_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 hubRtip::Float
: radius of the blade tipr::Array{Float}
: radial locations of interestaz::Float
: blade azimuth angleprecone::Float
: rotor precone angletilt::Float
: rotor tilt angle
FLOWFarm.get_peaks
— Methodget_peaks(array)
get the turning point values of a signal
Arguments
array::Array{Float}
: the signal to find the turning points, or peaks
FLOWFarm.get_peaks_indices
— Methodget_peaks_indices(array)
return the indices of the signal peaks
Arguments
array::Array{Float}
: the signal to find the turning points, or peaks
FLOWFarm.get_turb_atrbt_YAML
— Methodget_turb_atrbt_YAML(file_name)
read in turbine attributes from .yaml
Arguments
file_name::String
: path/to/attribute/file.yaml
FLOWFarm.get_turb_loc_YAML
— Methodget_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
FLOWFarm.get_wind_rose_YAML
— Methodget_wind_rose_YAML(file_name)
read in wind resource information from .yaml
Arguments
file_name::String
: path/to/wind/resource/file.yaml
FLOWFarm.grid_points
— Methodgrid_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
FLOWFarm.hermite_spline
— Methodhermite_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 yx0::Float
: x position of upwind endpoint of splinex1::Float
: x position of downwind endpoint of spliney0::Float
: y position of upwind endpoint of splinedy0::Float
: slope at upwind endpoint of spliney1::Float
: y position of downwind endpoint of splinedy1::Float
: slope at downwind endpoint of spline
FLOWFarm.iea37cs4BndryVRIntPM
— Methodiea37cs4BndryVRIntPM(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 boundarybndry_y_clsd::Array{Float,1}
: 1-D array of y-coordinates for the vertices around a singlar closed boundarybndry_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
FLOWFarm.latlong_to_xy
— Methodlatlong_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)
FLOWFarm.met2cart
— Methodmet2cart(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
FLOWFarm.multiple_components_op
— Functionmultiple_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 rV::Array{Float}
: v velocity component of the inflow at each rW::Array{Float}
: w velocity component of the inflow at each rOmega::Float
: rotor rotational speedr::Array{Float}
: radial locations of interestprecone::Float
: rotor precone angleyaw::Float
: rotor yaw angletilt::Float
: rotor tilt angleazimuth::Float
: blade azimuth anglerho::Float
: air density
Keyword Arguments
mu::Float
: air viscocity (can usually use the default)asound::Float
: speed of sound (can usually use the default)
FLOWFarm.nansafenorm
— Methodnansafenorm(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.
FLOWFarm.nansafesqrt
— Methodnansafesqrt(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()
FLOWFarm.point_velocity
— Methodpoint_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 interestturbine_x::Array{TF,nTurbines}
: turbine east-west locations in the state reference frameturbine_y::Array{TF,nTurbines}
: turbine north-south locations in the state reference frameturbine_z::Array{TF,nTurbines}
: turbine base height in the state reference frameturbine_yaw::Array{TF,nTurbines}
: turbine yaw for the given wind direction in radiansturbine_ct::Array{TF,nTurbines}
: turbine thrust coefficients for the given stateturbine_ai::Array{TF,nTurbines}
: turbine axial induction for the given staterotor_diameter::Array{TF,nTurbines}
: turbine rotor diametershub_height::Array{TF,nTurbines}
: turbine hub heightsturbine_local_ti::Array{TF,nTurbines}
: turbine local turbulence intensity for the given statesorted_turbine_index::Array{TF,nTurbines}
: array containing indices of wind turbines from most upwind to most downwind turbine in the given statewtvelocities::Array{TF,nTurbines}
: effective inflow wind speed for given statewind_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 1downwind_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)
FLOWFarm.pointinpolygon
— Functionpointinpolygon(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 interestvertices::Vector{Matrix{Number}(2)
: vertices of polygonnormals::Vector{Matrix{Number}(2)
: if not provided, they will be calculateds::Number
: smoothing factor for ksmax function (smoothmax)method::String
: currently only raycasting is availableshift::Float
: how far to shift point if it lies on an edge or vertexreturn_distance::Bool
: if true, return distance. if false, return -1 if in polygon or on the boundary, and 1 otherwise.
FLOWFarm.pointonline
— Methodpointonline(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 interestv1::Vector{Number}(2)
: first vertex of the linev2::Vector{Number}(2)
: second vertex of the linetol::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
FLOWFarm.rainflow
— Functionrainflow(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
FLOWFarm.ray_casting_boundary
— Methodray_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 faceturbine_x::Array{Float}
: turbine x locationsturbine_y::Array{Float}
: turbine y locationsdiscrete::Bool
: if true, indicates the boundary is made of multiple discrete regionss::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 discontinuityreturn_region::bool
: if true, return a vector specifying which region each turbine is in
FLOWFarm.recolor_jacobian!
— Methodrecolorjacobian!(sparsestruct::T,windstateid,nvariables,nturbines)
Helper function that recolors the jacobian for a single wind state
Arguments
sparse_struct
: sparseAEPstructunstablepattern structwind_state_id
: Wind state idn_variables
: Number of design variablesn_turbines
: Number of turbines
FLOWFarm.rediscretize_windrose
— Methodrediscretize_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 rosendirectionbins::Integer
: number of direction bins for the new wind rosestart::Float
: direction for first bin in radiansaveragespeed::Bool
: set whether or not to return the average wind speed as the speed for all bins
FLOWFarm.rotate_to_wind_direction
— Methodrotate_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 frameylocs::Array
: contains turbine north-south locations in the global reference framewind_direction_met::Array
: contains wind direction in radians in meteorological standard system (N=0 rad, proceeds CW, wind from direction given)
FLOWFarm.rotor_sample_points
— Functionrotor_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
: controls how many sample points to generatealpha::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 griduse_perimeter_points
: whether or not to include point exactly on the perimeter of the rotor swept area
FLOWFarm.round_farm_random_start
— Methodround_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 diametercenter::Number
: wind farm centerradius::Number
: wind farm radiusdiameter::Array{T,1}
: diameters of all wind turbines
FLOWFarm.smooth_max
— Methodsmooth_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 comparisony::Float
: second value for comparisons::Float
: controls the level of smoothing used in the smooth max
FLOWFarm.smooth_max
— Methodsmooth_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 valuess::Float
: controls the level of smoothing used in the smooth max
FLOWFarm.sparse_spacing!
— Methodsparsespacing!(spacingvec,turbinex,turbiney,relevant)
Helper function that calculates the relevant spacing constraints
Arguments
spacing_vec
: Vector containing the spacing constraintsturbine_x
: Vector containing x positions of turbinesturbine_y
: Vector containing y positions of turbinesrelevant
: 2d array that holds the relevant turbine pairs for the spacing constraint
FLOWFarm.splined_boundary!
— Methodsplined_boundary!(bndry_cons, 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
bndry_cons:: Vector{Float}
: boundary constraintsturbine_x::Array{Float}
: turbine x locationsturbine_y::Array{Float}
: turbine y locationsbndry_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"
FLOWFarm.splined_boundary
— Methodsplined_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 locationsturbine_y::Array{Float}
: turbine y locationsbndry_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"
FLOWFarm.sunflower_points
— Methodsunflower_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 generatealpha::Float
: Controls the smoothness of the boundary. alpha=0 is the standard "jagged edge" sunflower algoirthm and alpha=1 results in a smooth boundary.
FLOWFarm.turbine_powers_one_direction
— Methodturbine_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 onlyair_density::Float
power_models::Array{nturbines})
elements of array should be be of sub-types or AbstractPowerModel
FLOWFarm.turbine_spacing!
— Methodturbine_spacing!(spacing_vec,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
spacing_vec::Array{Float}
: vector of distances between turbinesturbine_x::Array{Float}
: turbine x locationsturbine_y::Array{Float}
: turbine y locations
FLOWFarm.turbine_spacing
— Methodturbine_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 locationsturbine_y::Array{Float}
: turbine y locations
FLOWFarm.turbine_velocities_one_direction
— Methodturbine_velocities_one_direction(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::Int=1, velocity_only::Bool=true, turbine_velocities=nothing,
turbine_ct=nothing, turbine_ai=nothing, turbine_local_ti=nothing)
Arguments
turbine_x::Array{TF,nTurbines}
: turbine east-west locations in the state reference frameturbine_y::Array{TF,nTurbines}
: turbine north-south locations in the state reference frameturbine_z::Array{TF,nTurbines}
: turbine base height in the state reference framerotor_diameter::Array{TF,nTurbines}
: turbine rotor diametershub_height::Array{TF,nTurbines}
: turbine hub heightsturbine_yaw::Array{TF,nTurbines}
: turbine yaw for the given wind direction in radianssorted_turbine_index::Array{TF,nTurbines}
: turbine sorted order upstream to downstream for given statect_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 analysiswind_farm_state_id::Int
: index to correct state to use from wind resource provided. Defaults to 1
FLOWFarm.unstable_sparse_aep_gradient!
— Methodunstablesparseaepgradient!(sparsestruct::T,x,farm,windstateid;prealloc_id=1,lock=nothing)
Function that calculates the AEP gradient for a single wind state using an unstable sparsity pattern
Arguments
sparse_struct
: sparseAEPstructunstablepattern structx
: Vector containing the design variablesfarm
: WindFarm structwind_state_id
: Wind state idprealloc_id
: Preallocation id (to select the correct preallocated memory inside the wind farm struct)lock
: SpinLock object to lock the farm struct for multithreadeding
FLOWFarm.update_safe_design_variables!
— Methodupdatesafedesignvariables!(spacingstruct::T,x)
Function that updates the safe design variables for the spacing constraints
Arguments
spacing_struct
: sparsespacingstruct structx
: Vector containing the design variables
FLOWFarm.update_turbine_powers!
— Methodupdateturbinepowers!(sparse_struct::T,i)
Helper function that updates the turbine powers for a single wind state from the sparse struct
Arguments
sparse_struct
: sparseAEPstructstablepattern structi
: Wind state id
FLOWFarm.wake_count_iec
— Methodwake_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 frameturbiney::Array{T,1}
: y locations of turbines in global reference framewinddirection::Float
: wind direction in radians in meteorological coordinates (0 rad. = from North)diameter::Array{T,1}
: diameters of all wind turbines
FLOWFarm.wake_deficit_model
— Methodwake_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, wake_deficits, contribution_matrix, deflections, current_index_loop, wind_speed_internal, sigma_squared, wtvelocities, sorted_turbine_index, model::CumulativeCurl)
Computes the wake deficit at a given location using the Cumulative Curl Model https://doi.org/10.5194/wes-2022-17
Arguments
locx::Float
: x coordinate where wind speed is calculatedlocy::Float
: y coordinate where wind speed is calculatedlocz::Float
: z coordinate where wind speed is calculatedturbine_x::Array(Float)
: vector containing x coordinates for all turbines in farmturbine_y::Array(Float)
: vector containing y coordinates for all turbines in farmturbine_z::Array(Float)
: vector containing z coordinates for all turbines in farmdeflection_y::Float
: deflection in the y direction of downstream wakedeflection_z::Float
: deflection in the z direction of downstream wakeupstream_turbine_id::Int
: index of the upstream wind turbine creating the wakedownstream_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 farmrotor_diameter::Array(Float)
: vector containing rotor diameters for all turbines in farmturbine_ai::Array(Float)
: vector containing initial velocity deficits for all turbines in farmturbine_local_ti::Array(Float)
: vector containing local turbulence intensities for all turbines in farmturbine_ct::Array(Float)
: vector containing thrust coefficients for all turbines in farmturbine_yaw::Array(Float)
: vector containing the yaw angle? for all turbines in farmwake_deficits
: matrix containing the wake deficits from every turbine to every other turbine for use in sparsity codescontribution_matrix
: matrix containing the contribution coefficients defined in the Cumlative Curl Modelmodel::CumulativeCurl
: indicates the wake model in use
FLOWFarm.wake_deficit_model
— Methodwake_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 calculatedlocy::Float
: y coordinate where wind speed is calculatedlocz::Float
: z coordinate where wind speed is calculatedturbine_x::Array(Float)
: vector containing x coordinates for all turbines in farmturbine_y::Array(Float)
: vector containing y coordinates for all turbines in farmturbine_z::Array(Float)
: vector containing z coordinates for all turbines in farmdeflection_y::Float
: deflection in the y direction of downstream wakedeflection_z::Float
: deflection in the z direction of downstream wakeupstream_turbine_id::Int
: index of the upstream wind turbine creating the wakedownstream_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 farmrotor_diameter::Array(Float)
: vector containing rotor diameters for all turbines in farmturbine_ai::Array(Float)
: vector containing initial velocity deficits for all turbines in farmturbine_local_ti::Array(Float)
: vector containing local turbulence intensities for all turbines in farmturbine_ct::Array(Float)
: vector containing thrust coefficients for all turbines in farmturbine_yaw::Array(Float)
: vector containing the yaw angle? for all turbines in farmmodel::GaussOriginal
: indicates the wake model in use
FLOWFarm.wake_deficit_model
— Methodwake_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, wake_deficits, contribution_matrix, deflections, current_index_loop, wind_speed_internal, sigma_squared, wtvelocities, sorted_turbine_index, 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 calculatedlocy::Float
: y coordinate where wind speed is calculatedlocz::Float
: z coordinate where wind speed is calculatedturbine_x::Array(Float)
: vector containing x coordinates for all turbines in farmturbine_y::Array(Float)
: vector containing y coordinates for all turbines in farmturbine_z::Array(Float)
: vector containing z coordinates for all turbines in farmdeflection_y::Float
: deflection in the y direction of downstream wakedeflection_z::Float
: deflection in the z direction of downstream wakeupstream_turbine_id::Int
: index of the upstream wind turbine creating the wakedownstream_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 farmrotor_diameter::Array(Float)
: vector containing rotor diameters for all turbines in farmturbine_ai::Array(Float)
: vector containing initial velocity deficits for all turbines in farmturbine_local_ti::Array(Float)
: vector containing local turbulence intensities for all turbines in farmturbine_ct::Array(Float)
: vector containing thrust coefficients for all turbines in farmturbine_yaw::Array(Float)
: vector containing the yaw angle? for all turbines in farmmodel::GaussSimple
: indicates the wake model in use
FLOWFarm.wake_deficit_model
— Methodwake_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, wake_deficits, contribution_matrix, deflections, current_index_loop, wind_speed_internal, sigma_squared, wtvelocities, sorted_turbine_index, 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 calculatedlocy::Float
: y coordinate where wind speed is calculatedlocz::Float
: z coordinate where wind speed is calculatedturbine_x::Array(Float)
: vector containing x coordinates for all turbines in farmturbine_y::Array(Float)
: vector containing y coordinates for all turbines in farmturbine_z::Array(Float)
: vector containing z coordinates for all turbines in farmdeflection_y::Float
: deflection in the y direction of downstream wakedeflection_z::Float
: deflection in the z direction of downstream wakeupstream_turbine_id::Int
: index of the upstream wind turbine creating the wakedownstream_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 farmrotor_diameter::Array(Float)
: vector containing rotor diameters for all turbines in farmturbine_ai::Array(Float)
: vector containing initial velocity deficits for all turbines in farmturbine_local_ti::Array(Float)
: vector containing local turbulence intensities for all turbines in farmturbine_ct::Array(Float)
: vector containing thrust coefficients for all turbines in farmturbine_yaw::Array(Float)
: vector containing the yaw angle? for all turbines in farmmodel::GaussYawVariableSpread
: indicates the wake model in use
FLOWFarm.wake_deficit_model
— Methodwake_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 calculatedlocy::Float
: y coordinate where wind speed is calculatedlocz::Float
: z coordinate where wind speed is calculatedturbine_x::Array(Float)
: vector containing x coordinates for all turbines in farmturbine_y::Array(Float)
: vector containing y coordinates for all turbines in farmturbine_z::Array(Float)
: vector containing z coordinates for all turbines in farmdeflection_y::Float
: deflection in the y direction of downstream wakedeflection_z::Float
: deflection in the z direction of downstream wakeupstream_turbine_id::Int
: index of the upstream wind turbine creating the wakedownstream_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 farmrotor_diameter::Array(Float)
: vector containing rotor diameters for all turbines in farmturbine_ai::Array(Float)
: vector containing initial velocity deficits for all turbines in farmturbine_local_ti::Array(Float)
: vector containing local turbulence intensities for all turbines in farmturbine_ct::Array(Float)
: vector containing thrust coefficients for all turbines in farmturbine_yaw::Array(Float)
: vector containing the yaw angle? for all turbines in farmmodel::GaussYaw
: indicates the wake model in use
FLOWFarm.wake_deficit_model
— Methodwake_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 calculatedlocy::Float
: y coordinate where wind speed is calculatedlocz::Float
: z coordinate where wind speed is calculatedturbine_x::Array(Float)
: vector containing x coordinates for all turbines in farmturbine_y::Array(Float)
: vector containing y coordinates for all turbines in farmturbine_z::Array(Float)
: vector containing z coordinates for all turbines in farmdeflection_y::Float
: deflection in the y direction of downstream wakedeflection_z::Float
: deflection in the z direction of downstream wakeupstream_turbine_id::Int
: index of the upstream wind turbine creating the wakedownstream_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 farmrotor_diameter::Array(Float)
: vector containing rotor diameters for all turbines in farmturbine_ai::Array(Float)
: vector containing initial velocity deficits for all turbines in farmturbine_local_ti::Array(Float)
: vector containing local turbulence intensities for all turbines in farmturbine_ct::Array(Float)
: vector containing thrust coefficients for all turbines in farmturbine_yaw::Array(Float)
: vector containing the yaw angle? for all turbines in farmmodel::JensenCosine
: indicates the wake model in use
FLOWFarm.wake_deficit_model
— Methodwake_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 calculatedlocy::Float
: y coordinate where wind speed is calculatedlocz::Float
: z coordinate where wind speed is calculatedturbine_x::Array(Float)
: vector containing x coordinates for all turbines in farmturbine_y::Array(Float)
: vector containing y coordinates for all turbines in farmturbine_z::Array(Float)
: vector containing z coordinates for all turbines in farmdeflection_y::Float
: deflection in the y direction of downstream wakedeflection_z::Float
: deflection in the z direction of downstream wakeupstream_turbine_id::Int
: index of the upstream wind turbine creating the wakedownstream_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 farmrotor_diameter::Array(Float)
: vector containing rotor diameters for all turbines in farmturbine_ai::Array(Float)
: vector containing initial velocity deficits for all turbines in farmturbine_local_ti::Array(Float)
: vector containing local turbulence intensities for all turbines in farmturbine_ct::Array(Float)
: vector containing thrust coefficients for all turbines in farmturbine_yaw::Array(Float)
: vector containing the yaw angle? for all turbines in farmmodel::JensenTopHat
: indicates the wake model in use
FLOWFarm.wake_deficit_model
— Methodwake_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 calculatedlocy::Float
: y coordinate where wind speed is calculatedlocz::Float
: z coordinate where wind speed is calculatedturbine_x::Array(Float)
: vector containing x coordinates for all turbines in farmturbine_y::Array(Float)
: vector containing y coordinates for all turbines in farmturbine_z::Array(Float)
: vector containing z coordinates for all turbines in farmdeflection_y::Float
: deflection in the y direction of downstream wakedeflection_z::Float
: deflection in the z direction of downstream wakeupstream_turbine_id::Int
: index of the upstream wind turbine creating the wakedownstream_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 farmrotor_diameter::Array(Float)
: vector containing rotor diameters for all turbines in farmturbine_ai::Array(Float)
: vector containing initial velocity deficits for all turbines in farmturbine_local_ti::Array(Float)
: vector containing local turbulence intensities for all turbines in farmturbine_ct::Array(Float)
: vector containing thrust coefficients for all turbines in farmturbine_yaw::Array(Float)
: vector containing the yaw angle? for all turbines in farmmodel::MultiZone
: indicates the wake model in use
FLOWFarm.wake_deflection_model
— Methodwake_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"
FLOWFarm.wake_deflection_model
— Methodwake_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"
FLOWFarm.wake_deflection_model
— Methodwake_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"
FLOWFarm.wake_deflection_model
— Methodwake_deflection_model(locx, locy, locz, turbine_id, turbine_definition::TurbineDefinition, model::NoYawDeflection, windfarmstate::SingleWindFarmState)
Bypasses yaw deflection calculations.
FLOWFarm.write_turb_loc_YAML
— Methodwrite_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