OWENSOpenFASTWrappers

Types and functions

OWENSOpenFASTWrappers.EnvironmentType

Environment(rho::TF1, dt::TF2, steplast::TAI) Environment(rho) = Environment(rho)

Contains specications for turbine environment/operating conditions as well as some backend memory

Inputs

  • rho::TF1: Working fluid density (kg/m^3)
  • dt::TF2: timestep for ADI
  • num_channels::TI1: number of output channels from AD15
  • steplast::TAI: prior simulation step index, used for unsteady wake propogation

Outputs:

  • none:
OWENSOpenFASTWrappers.StructureType

Structure

Contains the position and orientation info passed to AD15. These positions are all in global and include the turbine RefPos.

Inputs:

  • hubPos::TAF1:
  • hubOrient::TAF2:
  • hubVel::TAF3:
  • hubAcc::TAF4:
  • nacPos::TAF5:
  • nacOrient::TAF6:
  • nacVel::TAF7:
  • nacAcc::TAF8:
  • rootPos::TAF9:
  • rootOrient::TAF10:
  • rootVel::TAF11:
  • rootAcc::TAF12:
  • meshPos::TAF13:
  • meshOrient::TAF14:
  • meshVel::TAF15:
  • meshAcc::TAF16:

Outputs:

  • none:
OWENSOpenFASTWrappers.adiCalcOutputMethod
adiCalcOutput( )

calls AeroDynInflowC_CalcOutput to calculate resulting loads. Call this only after SetRotorMotion on all rotors/turbines.

Inputs:

time::c_double: current timestep

  • num_channels::int: number of output channels
OWENSOpenFASTWrappers.adiGetRotorLoadsMethod
adiGetRotorLoads(iturb; )

Gets the loads for a single turbine rotor

Inputs:

  • iturb::int: required, current turbine number

Outputs:

  • meshFrcMom::Array(float): loads from ADI at mesh nodes
OWENSOpenFASTWrappers.adiInitMethod
adiInit( output_root_name ; )

calls aerodyninflowinit to initialize AeroDyn and InflowWind together

Inputs:

  • ad_input_file_passed::int: flag to indicate the AD15 input file is passed as a string 0=false, 1=true (set to false if passing input file name instead, NOT SUPPORTED YET)

  • ad_input_file::string: name of input file for AD15 – this is read by julia and passed to AD15

  • ifw_input_file_passed::int: flag to indicate the InflowWind input file is passed as a string 0=false, 1=true (set to false if passing input file name instead, NOT SUPPORTED YET)

  • ifw_input_file::string: name of input file for InflowWind – this is read by julia and passed to InflowWind

  • gravity::float: optional, gravity value (default: 9.80665 m/s^2)

  • defFldDens::float: optional, air density (default: 1.225 kg/m^3)

  • defKinVisc::float: optional, kinematic viscosity of working fluid (default: 1.464E-05 m^2/s)

  • defSpdSound::float: optional, speed of sound in working fluid (default: 335.0 m/s)

  • defPatm::float: optional, atmospheric pressure (default: 103500.0 Pa) [used only for an MHK turbine cavitation check]

  • defPvap::float: optional, vapour pressure of working fluid (default: 1700.0 Pa) [used only for an MHK turbine cavitation check]

  • WtrDpth::float: optional, water depth (default: 0.0 m) [used only for an MHK turbine]

  • MSL2SWL::float: optional, offset between still-water level and mean sea level (default: 0.0 m) [positive upward, used only for an MHK turbine]

  • storeHHVel::int: optional, internal parameter for adi_library. Exposed for convenience, but not needed. [0=false, 1=true]

  • WrVTK::int: optional, write VTK output files at all timesteps to visualize AeroDyn 15 meshes [0 none (default), 1 ref, 2 motion]

  • WrVTK_Type::int: optional, write VTK output files as [1 surfaces (default), 2 lines, 3 both]

  • VTKNacDim::Array(float*6) optional, Nacelle Dimension for VTK visualization x0,y0,z0,Lx,Ly,Lz (m)

  • VTKHubRad::float: optional, HubRadius for VTK visualization (m)

  • wrOuts::int: optional, file format for writing outputs [0 none (default), 1 txt, 2 binary, 3 both]

  • DT_Outs::float64: optional, timestep for outputs to file [0.0 (default) for every timestep]

  • interp_order::int: optional, interpolation order used internally [1 first order (default), 2 second order]

  • dt::float64: required, timestep for AD15 (needed for setting internal constants)

  • t_max::float: required, total expected simulation time – used only for setting VTK counter width

Outputs:

  • num_channels::int: number of output channels
  • channel_names::string: string of output channel names from ADI
  • channel_units::string: string of output channel units from ADI
OWENSOpenFASTWrappers.adiPreInitMethod
adiPreInit(adilib_filename numTurbines transposeDCM)

Does some pre-initializing of the ADI library to setup arrays for each turbine

Inputs:

  • adilib_filename::string: path and name of AeroDyn-Inflow dynamic library
  • numTurbines::int: required, number of turbines to setup ADI to handle
  • transposeDCM::int: required, transpose DCM internally in ADI to match calling code convention for direction cosine matrices (default: 1==true)
OWENSOpenFASTWrappers.adiSetRotorMotionMethod
adiSetRotorMotion(iturb; )

Sets the motions for a single turbine rotor

Inputs:

  • iturb::int: required, current turbine number

  • HubPos::Array(float): required, (x,y,z) position of hub

  • HubOrient::Array(float): required, orientation of hub as 9 element vector of flattened DCM

  • HubVel::Array(float): required, (TVx,TVy,TVz,RVx,RVy,RVz) velocity of hub, does not include rotational velocity, so this is extra like from a platform

  • HubAcc::Array(float): required, (TAx,TAy,TAz,RAx,RAy,RAz) acceleration of hub, does not include rotational accel, so this is extra like from a platform

  • NacPos::Array(float): required, (x,y,z) position of nacelle

  • NacOrient::Array(float): required, orientation of nacelle as 9 element vector of flattened DCM

  • NacVel::Array(float): required, (TVx,TVy,TVz,RVx,RVy,RVz) velocity of nacelle

  • NacAcc::Array(float): required, (TAx,TAy,TAz,RAx,RAy,RAz) acceleration of nacelle

  • RootPos::Array(float): required, size (numBlades,3) position vectors of roots

  • RootOrient::Array(float): required, size (numBlades,9) orientation DCMs flattened to array of 9 element vectors

  • RootVel::Array(float): required, size (numBlades,6) velocity vectors of roots

  • RootAcc::Array(float): required, size (numBlades,6) acceleration vectors of roots

  • numMeshNodes::Array(int): required, number of structural mesh points (total across all blades)

  • MeshPos::Array(float): required, size (sum(numMeshNodes),3) position vectors of mesh points

  • MeshOrient::Array(float): required, size (sum(numMeshNodes),9) orientation DCMs flattened to array of 9 element vectors

  • MeshVel::Array(float): required, size (sum(numMeshNodes),6) velocity vectors of mesh points

  • MeshAcc::Array(float): required, size (sum(numMeshNodes),6) acceleration vectors of mesh points

Outputs:

  • none:
OWENSOpenFASTWrappers.adiSetupRotorMethod
adiSetupRotor(iturb; )

Sets the initial locations of a single rotor (root orientations/positions etc)

Inputs:

  • iturb::int: required, current turbine number

  • isHAWT::bool: required, false: VAWT or cross-flow turbine, true: HAWT

  • intTurbPos::Array(float): required, (x,y,z) position of turbine

  • initHubPos::Array(float): required, (x,y,z) position of hub

  • initHubOrient::Array(float): required, orientation of hub as 9 element vector of flattened DCM

  • initNacellePos::Array(float): required, (x,y,z) position of nacelle

  • initNacelleOrient::Array(float): required, orientation of nacelle as 9 element vector of flattened DCM

  • numBlades::int: required, number of blades

  • initRootPos::Array(float): required, size (numBlades,3) position vectors of roots

  • initRootOrient::Array(float): required, size (numBlades,9) orientation DCMs flattened to array of 9 element vectors

  • numMeshNodes::Array(int): required, number of structural mesh points (total across all blades)

  • initMeshPos::Array(float): required, size (sum(numMeshNodes),3) position vectors of mesh points

  • initMeshOrient::Array(float): required, size (sum(numMeshNodes),9) orientation DCMs flattened to array of 9 element vectors

OWENSOpenFASTWrappers.adiUpdateStatesMethod
adiUpdateStates( )

calls AeroDynInflowC_UpdateStates to step ADI forward to the next timestep Call this only after SetRotorMotion on all rotors/turbines.

Inputs:

time::c_double: current timestep

time_next::c_double: next timestep

OWENSOpenFASTWrappers.advanceAD15Method
advanceAD15(t_new;ts=2*pi/(turbine.omega[1]*turbine.ntheta))

Runs a previously initialized aero model (see ?setupTurb) in the unsteady mode (can be repeateadly called, or called for a specific time, or repeatedly called for sections of time)

Inputs

  • t_new::float: new time (s); will run from last time specified from the last call, to the current time specified, or from t=ts if the first time called
  • mesh::: OWENSFEA mesh for the turbine structure
  • azi::: hub azimuth (radians)
  • dt::float: optional timestep

Outputs:

  • n_steps: number timesteps taken
  • Fx: Array(sum(numMeshNodes),ntheta) Turbine Fx (N)
  • Fy: Array(sum(numMeshNodes),ntheta) Turbine Fy (N)
  • Fz: Array(sum(numMeshNodes),ntheta) Turbine Fz (N)
  • Mx: Array(sum(numMeshNodes),ntheta) Turbine Mx (N-m)
  • My: Array(sum(numMeshNodes),ntheta) Turbine My (N-m)
  • Mz: Array(sum(numMeshNodes),ntheta) Turbine Mz (N-m)
OWENSOpenFASTWrappers.createGeneralTransformationMatrixMethod

createGeneralTransformationMatrix(angleArray,axisArray)

Calculates the transformation matrix assocaited with a general Euler rotation sequence.

#Input

  • angleArray: = array of angles for Euler rotation sequence
  • axisArray: = array of axis of rotations for Euler rotation

#Output

  • dcmTotal: = transformation matrix of specified euler rotation sequence
OWENSOpenFASTWrappers.deformAD15Method
deformAD15(u_j,udot_j,uddot_j,azi,Omega_rad,OmegaDot_rad,hubPos,hubAngle,hubVel,hubAcc)

Sets the inputs for AD15.

Inputs

  • u_j: mesh displacements – in hub coordinates, (m,rad)
  • udot_j: mesh velocity – in hub coordinates, (m/s,rad/s)
  • uddot_j: mesh velocity – in hub coordinates, (m/s,rad/s)
  • azi: current azimuth (rad)
  • Omega_rad: angular velocity of hub about hub axis (rad/s)
  • OmegaDot_rad: angular acceleration of hub about hub axis (rad/s^2)
  • hubPos: current global hubPos (x,y,z) vector (m)
  • hubAngle: 3 angle set for hub orientation (rad), no rotation from spinning
  • hubVel: hub velocity in global coords, 6-vector (m/s,rad/s), does not include rotational velocity, so this is extra like from a platform
  • hubAcc: hub acceleration in global coords, 6-vector (m/s^2,rad/s^2), does not include rotational accel, so this is extra like from a platform
OWENSOpenFASTWrappers.frame_convertMethod

frameconvert(initframevals, transmat)

Internal, transfers 6 DOFs element-wise to a new reference frame

Input

  • init_frame_vals::Vector{<:float}: Values in 6 degrees of freedom in the initial reference frame
  • trans_mat::Array{<:float}: Transformation matrix to the output reference frame

Output

  • out_frame_vals: Values in 6 degrees of freedom in the output reference frame
OWENSOpenFASTWrappers.getAD15MeshDCMMethod

getAD15MeshDCM(turbine,u_j,azi,hubAngle)

Extract the mesh points orientations for all the AD15 blades ordering here is important 1. root to tip of blades, in blade order 2. root to tip of bottom struts, in blade order 3. root to tip of top struts, in blade order

Inputs

  • turbine: turbine data storage
  • u_j: mesh displacements – in hub coordinates, (m,rad)
  • azi: current azimuth (rad)
  • hubAngle: 3 angle set for hub orientation (rad), no rotation from spinning

#FIXME: add averaging of orientations to get nodes within blade/strut

OWENSOpenFASTWrappers.getAD15MeshPosMethod

getAD15MeshPos(turbine,u_j,azi,nacPos,hubPos,hubAngle)

Extract the mesh points for all the AD15 blades ordering here is important 1. root to tip of blades, in blade order 2. root to tip of bottom struts, in blade order 3. root to tip of top struts, in blade order

Inputs

  • turbine: turbine data storage
  • u_j: mesh displacements – in hub coordinates, (m,rad)
  • azi: current azimuth (rad)
  • hubPos: current global hubPos (x,y,z) vector (m)
  • nacPos: current global nacPos (x,y,z) vector (m)
  • hubAngle: 3 angle set for hub orientation (rad) , no rotation from spinning
OWENSOpenFASTWrappers.getAD15MeshVelAccMethod

getAD15MeshVelAcc(turbine,meshPos,udotj,uddotj,azi,Omegarad,OmegaDotrad,nacPos,hubPos,hubAngle,hubVel,hubAcc)

Extract the mesh velocities and accelerations for all the AD15 blades ordering here is important 1. root to tip of blades, in blade order 2. root to tip of bottom struts, in blade order 3. root to tip of top struts, in blade order

Inputs

  • turbine: turbine data storage
  • rootPos: root positions from call to getAD15MeshPos
  • u_j: mesh displacements – in hub coordinates, (m,rad)
  • udot_j: mesh velocity – in hub coordinates, (m/s,rad/s)
  • uddot_j: mesh velocity – in hub coordinates, (m/s,rad/s)
  • azi: current azimuth (rad)
  • Omega_rad: angular velocity of hub about hub axis (rad/s)
  • OmegaDot_rad: angular acceleration of hub about hub axis (rad/s^2)
  • hubPos: current global hubPos (x,y,z) vector (m)
  • nacPos: current global nacPos (x,y,z) vector (m)
  • hubAngle: 3 angle set for hub orientation (rad), no rotation from spinning
  • hubVel: hub velocity in global coords (m/s,rad/s), does not include rotational velocity, so this is extra like from a platform
  • hubAcc: hub acceleration in global coords (m/s^2,rad/s^2), does not include rotational accel, so this is extra like from a platform
OWENSOpenFASTWrappers.getRootDCMMethod

getRootDCM(turbine,u_j,azi,hubAngle)

Note on angles The OWENSFEA mesh in OWENS uses +X as the blade/strut long axis. In AeroDyn, the blade axis is +Z. So for transforming the blades from OWENSFEA mesh coordinates, first rotate by -90 degrees about Y, then do the 3,2,1 coordinate transforms with (Twist,Theta,Psi) = (Rz,Ry,Rx) = (Yaw,Pitch,Roll) Psid – rotation about Z axis – Yaw (Rz) – degrees Thetad – rotation about Y axis – Pitch (Ry) – degrees Twist_d – rotation about X axis – Roll (Rx) – degrees The rotation sequence is Roll –> Pitch –> Yaw. In rotation matrix form, it is R = RzRyRx (a [3,2,1] matrix order).

Inputs

  • turbine: turbine data storage
  • u_j: mesh displacements – in hub coordinates, (m,rad)
  • azi: current azimuth (rad)
  • hubAngle: 3 angle set for hub orientation (rad), no rotation from spinning

#FIXME: add averaging of orientations to get nodes within blade/strut

OWENSOpenFASTWrappers.getRootPosMethod

getRootPos(turbine,u_j,azi,nacPos,hubPos,hubAngle)

Extract the root positions for all ADI blades

Inputs

  • turbine: turbine data storage
  • u_j: mesh displacements – in hub coordinates, (m,rad)
  • azi: current azimuth (rad)
  • nacPos: current global nacPos (x,y,z) vector (m)
  • hubPos: current global hubPos (x,y,z) vector (m)
  • hubAngle: 3 angle set for hub orientation (rad), no rotation from spinning
OWENSOpenFASTWrappers.getRootVelAccMethod

getRootVelAcc(turbine,rootPos,udotj,uddotj,azi,Omegarad,OmegaDotrad,nacPos,hubPos,hubAngle,hubVel,hubAcc)

Extract the root velocities and accelerations for all ADI blades

Inputs

  • turbine: turbine data storage
  • rootPos: root positions from call to getRootPos
  • azi: current azimuth (rad)
  • Omega_rad: angular velocity of hub about hub axis (rad/s)
  • OmegaDot_rad: angular acceleration of hub about hub axis (rad/s^2)
  • nacPos: current global nacPos (x,y,z) vector (m)
  • hubPos: current global hubPos (x,y,z) vector (m)
  • hubAngle: 3 angle set for hub orientation (rad) , no rotation from spinning
  • hubVel: hub velocity in global coords, 6-vector (m/s,rad/s), does not include rotational velocity, so this is extra like from a platform
  • hubAcc: hub acceleration in global coords, 6-vector (m/s^2,rad/s^2), does not include rotational accel, so this is extra like from a platform
OWENSOpenFASTWrappers.ifwcalcoutputMethod
ifwcalcoutput(position,time)

calls inflow wind clacoutput

Inputs

  • position::Array(float): x, y, z sample position (m)
  • time::float: sample time (s)

Outputs:

  • velocities: x, y, z velocity at sample position
OWENSOpenFASTWrappers.ifwinitMethod
ifwinit(inflowlib_filename ;HWindSpeed=6.87,turbsim_filename="path/test.bts")

calls inflow wind init

Inputs

  • inflowlib_filename::string: path and name of inflow-wind dynamic library
  • HWindSpeed::float: optional, backup steady windspeed (m/s)
  • turbsim_filename::string: path and name of turbsim data e.g. "path/test.bts"

Outputs:

  • none:
OWENSOpenFASTWrappers.setupTurbMethod

setupTurb(adilib,adinputfile,ifwinputfile,adirootname,bldx,bldz,B; rho = 1.225, gravity = 9.80665, defKinVisc = 1.464E-05, defSpdSound = 335.0, defPatm = 103500.0, defPvap = 1700.0, WtrDpth = 0.0, MSL2SWL = 0.0, storeHHVel = 0, #false transposeDCM= 1, #true WrVTK = 2, WrVTKType = 3, VTKNacDim = [-.10 ,-.10 ,-.10 ,.2 ,.2 ,.2], VTKHubRad = 1.0, adinstrut = 2, adidt = 0.05, aditmax = 10, adiwrOuts = 0, adiDT_Outs = 0.0, isHAWT = false, # false: VAWT or cross-flow turbine, true: HAWT numTurbines = 1)

Initializes aerodynamic models and sets up backend persistent memory to simplify intermittent calling within coupled solver loops

Inputs

  • adi_lib: path to adi library (.so, .dylib, .dll)
  • ad_input_file: input file for aerodyn15
  • ifw_input_file: input file for inflow wind
  • adi_rootname: rootname for vtk outputs
  • bld_x: Blade x shape
  • bld_z: Blade z shape
  • B: Number of blades
  • rho: working fluid density (kg/m^3)
  • gravity: Gravitational acceleration (m/s^2)
  • defKinVisc: Kinematic viscosity of working fluid (m^2/s)
  • defSpdSound: Speed of sound in working fluid (m/s)
  • defPatm: Atmospheric pressure (Pa) [used only for an MHK turbine cavitation check]
  • defPvap: Vapour pressure of working fluid (Pa) [used only for an MHK turbine cavitation check]
  • WtrDpth: Water depth (m)
  • MSL2SWL: Offset between still-water level and mean sea level (m) [positive upward]
  • storeHHVel: unused here
  • transposeDCM: 0=false, 1=true transpose DCM internally for calculations
  • WrVTK: write VTK files from adi to directory adi-vtk [0 none, 1 ref, 2 motion]
  • WrVTK_Type: write VTK files from adi to directory adi-vtk [1 surface, 2 lines, 3 both]
  • VTKNacDim: Nacelle Dimension for VTK visualization x0,y0,z0,Lx,Ly,Lz (m)
  • VTKHubRad: HubRadius for VTK visualization (m)
  • adi_wrOuts: file format to write to
  • adi_DT_Outs: output timestep to write at
  • adi_nstrut: createmeshstruts is hard coded for 2 struts per blade
  • adi_dt: timestep
  • adi_tmax: maximum time
  • hubPos: hub position in global coordinates, 3-vector (m). NOTE: AD15 assumes a different hub location than OWENS
  • hubAngle: hub axis angle, 3-vector (deg), no rotation from spinning
  • nacPos: nacelle position in global coordinates, 3-vector (m). NOTE: AD15 assumes a different hub location than OWENS
  • nacAngle: nacelle axis angle, 3-vector (deg)
  • numTurbines: number of turbines
  • isHAWT: # false: VAWT or cross-flow turbine, true: HAWT

Outputs:

  • none:
OWENSOpenFASTWrappers.transMatMethod

transMat(theta1, theta2, theta3)

Internal, computes the 3x3 transformation matrix for given input rotations. The generated matrix is the closest orthonormal matrix to the Bernoulli-Euler transformation matrix from beam theory, which assumes small rotations. A full description of this matrix is found in the "FASTCoordinateSystems.doc" document by Jason Jonkman.