OWENSFEA

Types and functions

OWENSFEA.ConcNDL1DType

Internal, NodalTerms node number, local dof (diagonal), and value

OWENSFEA.ConcNDL2DType

Internal, NodalTerms node number, local dof1, local dof2, and value

OWENSFEA.DispDataType

Internal, displacement, velocity, and acceleration for each element

OWENSFEA.DispOutType
DispOut(elStrain,displ_sp1,displddot_sp1,displdot_sp1)

Internal, displacement, velocity, and acceleration for each element

Inputs

  • elStrain: Not used, should be removed from this struct
  • `displ_sp1::Array{<:float}: displacement position for each dof
  • `displddot_sp1::Array{<:float}: displacement acceleration for each dof
  • `displdot_sp1::Array{<:float}: displacement velocity for each dof
OWENSFEA.ElType

Internal, see ?Ort and ?SectionPropsArray

OWENSFEA.ElInputType

Internal, inputs to Timoshenko element, inputs given by FEAModel struct and mesh

OWENSFEA.ElStrainType
ElStrain(eps_xx_0,eps_xx_z,eps_xx_y,gam_xz_0,gam_xz_y,gam_xy_0,gam_xy_z)

Struct containing element straing

Inputs

  • epsilon_x::float: epsilon_x strain in the x direction
  • epsilon_y::float: epsilon_y strain in the y direction
  • epsilon_z::float: epsilon_z strain in the z direction
  • kappa_x::float: kappa_x curvature in the x direction
  • kappa_y::float: kappa_y curvature in the y direction
  • kappa_z::float: kappa_z curvature in the z direction

Outputs:

  • none:
OWENSFEA.FEAModelMethod
FEAModel(;analysisType = "TNB",
    initCond = [],
    aeroElasticOn = false,
    guessFreq = 0.0,
    airDensity=1.2041,
    gravityOn = true,
    nlOn = false,
    spinUpOn = false,
    outFilename = "none",
    RayleighAlpha = 0.0,
    RayleighBeta = 0.0,
    elementOrder = 1,
    joint = [0,0],
    platformTurbineConnectionNodeNumber = 1,
    jointTransform = 0.0,
    reducedDOFList = zeros(Int,2),
    numDOFPerNode = 6,
    numNodes = 0,
    numModes = 20,
    nlParams = 0,
    pBC = 0,
    nodalTerms = 0.0,
    iterationType = "NR",
    adaptiveLoadSteppingFlag = true,
    tolerance = 1.0000e-06,
    maxIterations = 50,
    maxNumLoadSteps = 20,
    minLoadStepDelta = 0.0500,
    minLoadStep = 0.0500,
    prescribedLoadStep = 0.0)

Model inputs for FEA analysis, struct

Inputs

  • analysisType::string: Newmark Beta time stepping "TNB", Dean time stepping "TD", modal "M", and stiff "stiff" (where the forces are just directly mapped, the displacements and strains set to 0 and the structures not run)
  • initCond::Array{<:float}: Initial conditions Nx3 matrix consisting of nodeNumber, local DOF (1-6), and displacement value
  • aeroElasticOn::Bool: Include simplified flutter calculataions in the timoshenko element matrices
  • guessFreq::float: aeroelastic starting guess, only used if aeroElasticOn
  • airDensity::float: working fluid density
  • gravityOn::Bool orArray{<:float}: vector of 3 or flag to include distributed gravity acceleration (9.81m/s) in the negative z-direction
  • nlOn::Bool: flag for solver to calculate deflection induced stiffness changes and associated convergance to the coupled solution
  • spinUpOn::Bool: flag to perform a static analysis (warm start) prior to performing modal analysis
  • outFilename::string: /path/to/desired/output/filename if it doesn't exist already it is created, if exists, is overwritten
  • RayleighAlpha::float: Rayleigh alpha damping used in timoshenko beam damping matrix
  • RayleighBeta::float: Rayleigh beta damping used in timoshenko beam damping matrix
  • elementOrder::int: order of element: 1 linear, 2 quadratic
  • joint::Array{<:float}: jointNumber masterNode slaveNode jointType (0 weld/fixed, 1 pinned, 2 hinge along local "2", 3 hinge about local "1", 4 hinge along "3", 5 rigid bar constraint) jointMass 0.0 jointPsi jointTheta
  • platformTurbineConnectionNodeNumber::int: node at which reaction forces are calculated
  • jointTransform: not used as an input, is calculated, local transform between dependent and active DOFs for nodes associated with a particular joint
  • reducedDOFList::Array{<:int}: not used as an input, is calculated, map of original DOF numbering to reduced DOF numbering
  • numDOFPerNode::int: number of degrees of freedom per node
  • numNodes::int: total number of nodes in the mesh
  • numModes::int: number of modes to calculate
  • nlParams::NlParams: optional there in case the Nlparams struct is passed in, should be cleaned up since redundant
  • alpha::float64: optional newmark beta alpha parameter,If TD, use 0.25
  • gamma::float64: optional newmark beta gamma parameter, if static, use 0. If hydro, use 1.0
  • pBC::Array{<:float}: Nx3 array consisting of node, local dof, specified displacement value for the boundary condition
  • nodalTerms: Concentrated nodal terms, should be replaced with the nodal input data array and the calc done internally
  • iterationType::string: FEA displacement update calculation, Newton Raphson "NR", Direct Iteration "DI"
  • adaptiveLoadSteppingFlag: Unused, should be removed
  • tolerance::float: FEA total mesh unsteady analysis convergence tolerance for a timestep
  • maxIterations: FEA total mesh unsteady analysis convergence max iterations for a timestep
  • maxNumLoadSteps: used in static (steady state) analysis
  • minLoadStepDelta: used in static (steady state) analysis
  • minLoadStep: used in static (steady state) analysis
  • prescribedLoadStep: used in static (steady state) analysis
  • predef::Bool: will update the elStorage array if passed into Unsteady() with the nonlinear strain stiffening, to be used for subsequent analyses

Outputs:

  • none:
OWENSFEA.MeshType
Mesh(nodeNum,numEl,numNodes,x,y,z,elNum,conn,type,meshSeg,structuralSpanLocNorm,structuralNodeNumbers,structuralElNumbers)

Struct with mesh definition

Inputs

  • nodeNum::Array{<:int}: Number mapping of nodes (typically 1:Nnodes)
  • numEl::int: total number of elements
  • numNodes::int: total number of nodes
  • x::Array{<:float}: Nodal x position
  • y::Array{<:float}: Nodal y position
  • z::Array{<:float}: Nodal z position
  • elNum::Array{<:int}: Number mapping of elements (typically 1:Nelements)
  • conn::Array{<:int}: Nelemx2 connectivity between nodes, gaps between joints (which are defined in the joints)
  • type::Array{<:int}: 0-blade 1-tower 2-strut
  • meshSeg::Array{<:int}: number of nodes within each segment, with segments consisting of tower, blade 1 2 etc, struts
  • structuralSpanLocNorm::Array{<:float}: Should be named heigh loc norm - unitized position along the blade height, used for aeroload mapping
  • structuralNodeNumbers::Array{<:int}: Node numbers associated with blades for aero loads mapping
  • structuralElNumbers::Array{<:int}: Element numbers associated with blades for aero loads mapping
  • nonRotating::Array{<:int}: size(Nsections,numNodes) where nsections are the number of sections of the mesh that are non-rotating, like if for some reason you had two towers, or if you had multiple guy wires
  • hubNodeNum::int: Node number where the rotating part of the turbine starts, assumes meshing always starts with tower, then blades, etc.

Outputs:

  • none:
OWENSFEA.NlParamsType
NlParams(iterationType,adaptiveLoadSteppingFlag,tolerance,maxIterations,maxNumLoadSteps,minLoadStepDelta,minLoadStep,prescribedLoadStep)

See ?FEAModel

OWENSFEA.OrtType
Ort(Psi_d,Theta_d,Twist_d,Length,elNum,Offset)

Struct with element orientation

Inputs

  • Psi_d::Array{<:float}: length NumEl, element rotation about 3 in global FOR (deg) These angles are used to transform from the global coordinate frame to the local element/joint frame via a 3-2 Euler rotation sequence.
  • Theta_d::Array{<:float}: length NumEl, element rotation about 2 (deg)
  • Twist_d::Array{<:float}: length NumEl, element twist (deg)
  • Length::Array{<:float}: length NumEl, element length (m)
  • elNum::Array{<:float}: Element number the other arrays are associated with
  • Offset::Array{<:float}: hub frame coordinate of node 1 of the element

Outputs:

  • none:
OWENSFEA.SectionPropsArrayType
SectionPropsArray(ac,twist,rhoA,EIyy,EIzz,GJ,EA,rhoIyy,rhoIzz,rhoJ,zcm,ycm,a,EIyz,alpha1,alpha2,alpha3,alpha4,alpha5,alpha6,rhoIyz,b,a0,aeroCenterOffset,xaf,yaf)

Struct with element sectional properties, each component is a 1x2 array with distributed properties

Inputs

  • ac::Array{<:float}: aerodynamic center, used in flutter approximation
  • twist::Array{<:float}: element twist (rad)
  • rhoA::Array{<:float}: rho * A in standard SI units
  • EIyy::Array{<:float}: E * Iyy
  • EIzz::Array{<:float}: E * Izz
  • GJ::Array{<:float}: G * J
  • EA::Array{<:float}: E * A
  • rhoIyy::Array{<:float}: rho * Iyy
  • rhoIzz::Array{<:float}: rho * Izz
  • rhoJ::Array{<:float}: rho * J
  • zcm::Array{<:float}: z location of center of mass
  • ycm::Array{<:float}: y location of center of mass
  • a::Array{<:float}: possibly lift slope
  • EIyz::Array{<:float}: E * Iyz
  • alpha1::Array{<:float}: #This is always 0 in the element file, and it is unclear what it is used for since I can't find it being used in the code
  • alpha2::Array{<:float}: doesn't appear to be used
  • alpha3::Array{<:float}: doesn't appear to be used
  • alpha4::Array{<:float}: doesn't appear to be used
  • alpha5::Array{<:float}: doesn't appear to be used
  • alpha6::Array{<:float}: doesn't appear to be used
  • rhoIyz::Array{<:float}: rho * Iyz
  • b::Array{<:float}: used in flutter approximation, possibly a chord or thickness value
  • a0::Array{<:float}: zero lift angle of attack, used in flutter approximation
  • aeroCenterOffset::Array{<:float}: doesn't appear to be used
  • xaf::Array{<:float}: x airfoil coordinates (to scale)
  • yaf::Array{<:float}: y airfoil coordinates (to scale)

Outputs:

  • none:
OWENSFEA.ConcMassAssociatedWithElementMethod
ConcMassAssociatedWithElement(conn,joint,nodalMassTerms,nodalStiffnessTerms,nodalLoads)

Compiles concentrated mass, stiffness, and load associated with a node from both ndl and joint files. The mod* variables are passed back with these terms removed to prevent duplicate application of shared nodal terms between elements

#Input

  • conn connectivity list for element
  • joint joint array for nodal terms
  • nodalMassTerms listing of concentrated nodal mass terms
  • nodalStiffnessTerms listing of concentrated nodal stiffness terms
  • nodalLoads listing of concentrated nodal loads terms

#Output

  • mass array of concentrated mass associated with element
  • stiff array of concentrated stiffness associated with element
  • load array of concentrated loads associated with element
  • modJoint modified joint object removing nodal terms that have/will be applied to the element calculations
  • modNodalMassTerms modified nodal mass object removing nodal terms that have/will be applied to the element calculations
  • modalStiffnessTerms modified nodal stiffness object removing nodal terms that have/will be applied to the element calculations
  • modNodalLoads modified nodal loads object removing nodal terms that have/will be applied to the element calculations
OWENSFEA.ModalOutputMethod
writeOutput(freq,damp,phase1,phase2,imagComponentSign,fid)

Internal, writes an output file and or formats an output for modal analysis.

#Input

  • freq: array of modal frequencies
  • damp: array of modal damping ratios
  • phase1: array of in phase mode shapes
  • phase2: array of out of phase mode shapes
  • imagComponentSign: array of sign of imaginary components
  • fid: file identifier for output

#Output

  • freqSorted: array of sorted(by frequency) modal frequencies
  • dampSorted: array of sorted(by frequency) modal damping ratios
  • imagCompSignSorted: array of sorted(by frequency) of imaginarycomponentSign array
  • U_x_0: see ?Modal outputs
  • U_y_0: see ?Modal outputs
  • U_z_0: see ?Modal outputs
  • theta_x_0: see ?Modal outputs
  • theta_y_0: see ?Modal outputs
  • theta_z_0: see ?Modal outputs
  • U_x_90: see ?Modal outputs
  • U_y_90: see ?Modal outputs
  • U_z_90: see ?Modal outputs
  • theta_x_90: see ?Modal outputs
  • theta_y_90: see ?Modal outputs
  • theta_z_90: see ?Modal outputs
OWENSFEA.applyBCMethod
applyBC(Kg,Fg,BC,u,iterationType,numDofPerNode)

Internal, applies boundary conditions to the stiffness matrix and load vector for a static analysis.

#Input

  • Kg assembled global stiffness matrix
  • Fg assembled global load vector
  • BC struct of boundary condition information
  • u global displacement vector
  • iterationType for nonlinear analysis, not used in BLAST
  • numDofPerNode number of degrees of freedom per node

#Output

  • Kg global stiffness matrix with boundary conditions
  • Fg global load vector with boundary condition
OWENSFEA.applyBCModalMethod
applyBCModal(K,BC,numDofPerNode)

Internal, applies boundary conditions to a system matrix for modal analysis

Inputs

  • K: assembled global system matrix
  • BC: struct of boundary condition information
  • numDofPerNode: number of degrees of freedom per node

Outputs:

  • Knew global system matrix with boundary conditions
OWENSFEA.applyConcentratedTermsMethod
applyConcentratedTerms(numNodes, numDOFPerNode; filename="none",data=[1 "M6" 1 1 0.0], jointData=[])

Internal, applies 6x6 concentrated nodal terms from user input.

#Input

  • filename: string containing nodal terms filename
  • data: Nx5 or Nx4 array matching general [1 "M6" 1 1 0.0] or diagonal only [1 "M" 1 0.0] aligning with node, Type, dof, value, where type is M,C,K, or F

#Output

  • nodalTerms::NodalTerms: see ?NodalTerms object containing concentrated nodal data
OWENSFEA.applyConstraintsMethod

Internal, this function transforms a matrix by the transformation matrix to enforce joint constraints

OWENSFEA.applyConstraintsVecMethod

Internal, this function transforms a vector by the transformation matrix to enforce joint constraints

OWENSFEA.assembly!Method
assembly(Ke,Fe,conn,numNodesPerEl,numDOFPerNode,Kg,Fg)

Internal, assembles the element matrix and load vector into the global system of equations

#Input

  • Ke: element matrix
  • Fe: element vector
  • conn: element connectivity
  • numNodesPerEl: number of nodes per element
  • numDofPerNode: number of degrees of freedom per node
  • Kg: global system matrix
  • Fg: global load vector

#Output

  • Kg: global system matrix with assembled element
  • Fg: global load vector with assembled element
OWENSFEA.assemblyMatrixOnlyMethod
assemblyMatrixOnly(Ke,conn,numNodesPerEl,numDOFPerNode,Kg)

Internal, assembles the element matrix into the global system of equations

Inputs

  • Ke: element matrix
  • conn: element connectivity
  • numNodesPerEl: number of nodes per element
  • numDofPerNode: number of degrees of freedom per node
  • Kg: global system matrix

Outputs:

  • Kg: global system matrix with assembled element
OWENSFEA.autoCampbellDiagramMethod

frequencies = autoCampbellDiagram(FEAinputs,mymesh,myel,system,assembly,sections; minRPM = 0.0, maxRPM = 40.0, NRPM = 9, # int vtksavename = nothing, saveModes = [1,3,5], #must be int saveRPM = [1,3,5], #must be int mode_scaling = 500.0, )

Automated Campbell Diagram Generator, this function runs the model with centrifugal stiffening for the specified RPM levels. If FEAinputs.analysisType == "GX" and vtksavename are specified, it will output paraview mode shape files at the specified save name.

#Inputs

  • FEAinputs::OWENSFEA.FEAModel: The FEA modeling options
  • mymesh::OWENSFEA.Mesh: a previously generated turbine mesh
  • myel::OWENSFEA.El: the element properties associated with that mesh
  • system::GXBeam.System: the converted GXBeam system from the mesh and el
  • assembly::GXBeam.AssemblyState: the converted GXBeam assembly from the mesh and el
  • sections::Array{Float64, 3}: the 3D point cloud to be converted to VTK format
  • minRPM::Float64: minimum RPM to be run, e.x. 0.0
  • maxRPM::Float64: maximum RPM to be run e.x. 40.0
  • NRPM::Int64: number of linear discretizations of RPM e.x. 9 must be int
  • vtksavename::string: filename (with path if desired) of the VTK outputs if GX. Set to "nothing" to not save.
  • saveModes::Array{Int64}: The modes to save in the VTK outputs e.x. [1,3,5] must be int
  • saveRPM::Array{Int64}: The RPMs to save in the VTK outputs e.x. [1,3,5] must be int
  • mode_scaling::Float64: The mode scaling in the VTK outputs e.x. 500.0

#Outputs

  • frequency::Array{Float64}: The output modal frequencies
OWENSFEA.calcUnormMethod

This function calculates a relative norm between two vectors: unew and uold

OWENSFEA.calculateBCMapMethod

calculateBCMap(numpBC,pBC,numDofPerNode,reducedDofList)

Internal, creates a boundary condition map between full and reduced dof listing as a result of constraints.

#Input

  • numpBC number of boundary conditions
  • pBC array of boundary condition data
  • numDofPerNode number of degrees of freedom per node
  • reducedDofList array of reduced DOF numbering

#Output

  • elStorage map for boundary conditions between full and reduced dof list
OWENSFEA.calculateLambdaMethod
calculateLambda(theta1,theta2,theta3)

This function calculates a transformation matrix to transform the element degree of freedom vector (12 DOFs) from the hub frame to the element frame. The transformation matrix is constructed via the direction cosine matrices of a 3-2-1 Euler rotation sequence.

#Input *theta1::float: angle (rad) of rotation for 1st rotation of 3-2-1 sequence *theta2::float: angle (rad) of rotation for 2nd rotation of 3-2-1 sequence *theta3::float: angle (rad) of rotation for 3rd rotation of 3-2-1 sequence

#Output *lambda::Array{<:float}: 12 x 12 transformation matrix

OWENSFEA.calculateLoadVecFromDistForceMethod

calculateLoadVecFromDistForce(elementOrder,x,xloc,twist,sweepAngle,coneAngle,rollAngle,extDistF2Node,extDistF3Node,extDistF4Node)

Takes in a global 6dof distributed force at two nodal points and returns the 6dof force in the element FOR

#Input

  • elementOrder:::
  • x::Array{<:float}: mesh x-position
  • xloc::Array{<:float}: local x-position [0 elLength]
  • twist::Array{<:float}: element twist angle (rad)
  • sweepAngle::Array{<:float}: element sweep angle (rad)
  • coneAngle::Array{<:float}: element cone angle (rad)
  • rollAngle::Array{<:float}: element roll angle (rad)
  • extDistF2Node::Array{<:float}: turbine Tangential force
  • extDistF3Node::Array{<:float}: turbine Normal force
  • extDistF4Node::Array{<:float}: turbine M25 moment

#Output

  • Fe::Array{float}: 6x1 Force on element in element FOR
OWENSFEA.calculateROMMethod
calculateROM(model,mesh,el,displ,omegaVec,omegaDotVec,elStorage,countedNodes)

This function calculates a reduced order model for a conventional structural dynamics system (parked, non-rotating)

#Input

  • model object containing model data
  • mesh object containing mesh data
  • el object containing elementdata
  • displ displacement vector
  • rbData: vector containing rigid body displacement, velocity, and acceleration
  • elStorage object containing stored element data
  • countedNodes prevents applied nodal terms from double counting

#Output

  • rom object containing reduced order model data
OWENSFEA.calculateROMGyricMethod

calculateROMGyric(feamodel,mesh,el,displ,omegaVec,omegaDotVec,elStorage,rom0,countedNodes)

Calculates a reduced order feamodel with rotational/ rigid body motion effects

#Input

  • feamodel: object containing feamodel data
  • mesh: object containing mesh data
  • el: object containing elementdata
  • displ: displacement vector
  • rbData: vector of hub frame accel (1-3), angular velocity components (4-6), and angular accleration (7-9)
  • elStorage: object containing stored element data
  • rom0: object containing parked/conventional reduced order feamodel
  • countedNodes: prevents applied nodal terms from double counting

#Output

  • rom: object containing reduced order feamodel data
OWENSFEA.calculateReactionForceAtNodeMethod
calculateReactionForceAtNode(nodeNum,model,mesh,el,elStorage,timeInt,dispData,displ_iter,rbData,Omega,OmegaDot,CN2H,countedNodes)

Internal, calculates the reaction force at a node by post processing all element associated with a node through connectivity or joint constraints.

#Input

  • nodeNum: node number joint constraints are desired at
  • model: object containing model data
  • mesh: object containing mesh data
  • elStorage: object containing stored element data
  • el: object containing element data
  • timeInt: object containing time integration parameters
  • dispData: object containing displacement data
  • displ_iter: converged displacement solution
  • rbData: vector containing rigid body displacement, velocity, and acceleration
  • Omega: rotor speed (Hz)
  • OmegaDot: rotor acceleratin (Hz)
  • CN2H: transformation matrix from inertial frame to hub frame
  • countedNodes: prevents nodal terms from being double counted

#Output

  • cummulativeForce: vector containing reaction force at nodeNum
OWENSFEA.calculateReducedDOFVectorMethod

Internal, searches over all DOFs in a structural model and determines and returns "dofVector" containing only unconstrained DOFs

OWENSFEA.calculateShapeFunctionsMethod

calculateShapeFunctions(elementOrder,xi,x)

This function calculates the Lagrange shape function, shape function derivative, and Jacobian to map between the local element domain and physical length of the element. The shape function derivative is defined with respect to the physical length domain. The shape functions may be linear or quadratic in order.

#Input

  • elementOrder order of element: 1 linear, 2 quadratic
  • xi guass point values to evaluate shape functions at
  • x nodal coordinates in physical length domain

#Output

  • N shape function value at specified gauss points
  • p_N_x shape function derivative w.r.t physical length domain at specified gauss points
  • Jac Jacobian for mat between local element domain and physical length domain.
OWENSFEA.calculateStructureMassPropsMethod

calculateStructureMassProps(elStorage)

This function caclulates structural mass properties of the finite element mesh (mass, moment of inertia, mass center) about the origin of the mesh coordinate system.

#Input

  • elStorage::ElStorage see ?ElStorage, object containing arrays of stored element information

#Output

  • structureMass::float mass of structure
  • structureMOI::float moment of inertia tensor of structgure
  • structureMassCenter::float center of mass of structure
OWENSFEA.calculateTimoshenkoElementInitialRunMethod
calculateTimoshenkoElementInitialRun(elementOrder,modalFlag,xloc,sectionProps,sweepAngle,coneAngle,rollAngle,aeroSweepAngle,x,y,z,concMassFlag,concMass,Omega)

Internal, see ?initialElementCalculations, performs initial element calculations and stores them for later use and efficiency gains.

OWENSFEA.calculateTimoshenkoElementNLMethod
calculateTimoshenkoElementNL(input,elStorage;predef=nothing)

Internal, performs nonlinear element calculations.

#Inputs

  • input::ElInput: see ?ElInput
  • elStorage::ElStorage: see ?ElStorage
  • predef::Bool: optional, if true, mutates ElStorage to include the nonlinear strain stiffening

#Outputs

  • ElOutput: see ?ElOutput
OWENSFEA.calculateTimoshenkoElementNLSSMethod
calculateTimoshenkoElementNLSS(elinput)

Performs selective nonlinear element calculations. Only stiffness matrix contributions are evaluate. No other calculations are performed to facilitate efficiency gains.

#Input

  • elinput: object containing element input

#Output

  • eloutput: object containing element data
OWENSFEA.calculateTimoshenkoElementStrainMethod
calculateTimoshenkoElementStrain(elementOrder,nlOn,xloc,sectionProps,sweepAngle,coneAngle,rollAngle,aeroSweepAngle,disp)

Internal, calculates element strain for a Timoshenko element

#Outputs

  • ElStrain: See ?ElStrain
OWENSFEA.constructReducedDispVecFromEigVecMethod
constructReducedDispVecFromEigVec(vec1,reducedDOFList,BC)

Internal, This function takes the original mode shape and modifies it to account for boundary conditions

Inputs

  • vec1:
  • reducedDOFList:
  • BC:

Outputs:

  • vec1Red:
OWENSFEA.createJointTransformMethod

createJointTransform(joint,numNodes,numDofPerNode)

Internal, calculates the JointTransform of a structural system.

#Input

  • joint: object containing joint data
  • numModes: number of nodes in mesh
  • numDofPerNode: number of degrees of freedom per node

#Output

  • jointTransform: joint transformation matrix
  • reducedDOF: map of original DOF numbering to reduced DOF numbering
OWENSFEA.createTdaMethod

Internal, creates a constraint transformation matrix for a single joint. Tda is this matrix, dDOF contains a listing of dependent global DOFs associated with this joint, and aDOF contains a listing of active global DOFs associated with this joint.

OWENSFEA.extractFreqDampMethod
extractFreqDamp(val,vec,numDOFPerNode,jointTransform,reducedDOFList,BC,analysisType)

Internal, calculates the eigenvalues and vectors of a structural dynamic system

Inputs

  • val: eigenvalue
  • vec: eigenvector
  • numDOFPerNode: number of degrees of freedom per node
  • jointTransform: joint transformation matrix from reduced to full DOF list
  • reducedDOFList: listing of reduced DOFs
  • BC: boundary condition object containing boundary condition info
  • analysisType: analysis type

Outputs:

  • freq: modal frequency
  • damp: modal damping
  • phase1: in phase mode shape (real part of mode shape)
  • phase2: out of phase mode shape (imaginary part of mode shape)
  • sortedModes: total, complex mode shape
OWENSFEA.extractdaInfoMethod

Internal, gets the total number of DOFs in the model, active number of DOFs in the model, and a list of slave DOFs that will be eliminated by joint constraints.

OWENSFEA.findElementsAssociatedWithNodeNumberMethod
findElementsAssociatedWithNodeNumber(nodeNum,conn,jointData)

Internal, finds elements associated with a node number through mesh connectivity or joint constraints

#Input

  • nodeNum node number joint constraints are desired at
  • conn object containing mesh connectivity
  • jointData object containing joint information

#Output

  • elList array containing a list of element numbers associated with nodeNum
  • localNode array containing the local node number that correspond to nodeNum in the list of associated elements
OWENSFEA.getGPMethod
getGP(numGP)

Internal, defines gauss point coordinates in a local element frame and the associated weights for Gaussian quadrature numerical integration.

#Input

  • numGP: number of quad points used for integration

#Output

  • xi: list of quad point coordinates in local element frame
  • weight: associated weights for quad point coordinate
OWENSFEA.initialElementCalculationsMethod
initialElementCalculations(feamodel,el,mesh)

performs initial element calculation for use later in analysis for efficiency gains.

Inputs

  • feamodel::FEAmodel: see ?Feamodel
  • el::El: see ?El
  • mesh::Mesh: see ?Mesh

Outputs:

  • elStorage:ElStorage: see ?ElStorage
OWENSFEA.makeBCdataMethod
makeBCdata(pBC,numNodes,numDofPerNode,reducedDOFList,jointTransform)

Internal, usese the pBC matrix and calculates/stores boundary condition data

#Input

  • pBC See ?FEAModel.pBC
  • numNodes number of nodes in structural model
  • numDofPerNode number of degrees of freedom per node
  • reducedDOFList joint transformation matrix from reduced to full DOF list
  • jointTransform listing of reduced DOFs

#Output

  • BC:BC_struct see ?BC_struct
OWENSFEA.mapMatrixNonSymMethod

Internal, function to form total stifness matrix and transform to desired DOF mapping

OWENSFEA.mapVectorMethod

Internal, forms total force vector and transform to desired DOF mapping

OWENSFEA.modalMethod
modal(feamodel,mesh,el;Omega=0.0,displ=zeros(mesh.numNodes*6),OmegaStart=0.0,returnDynMatrices=false)

Modal analysis

Inputs

  • feamodel::FEAModel: see ?FEAModel
  • mesh::Mesh: see ?Mesh
  • el::El: see ?El
  • Omega::float: Rotational rate in Hz
  • displ::Array{<:float}: zeros(mesh.numNodes*6) initial (warm start) displacements for each dof
  • OmegaStart::float: rotor speed (Hz) from previous analysis if stepping through various rotor speeds, may be useful in load stepping
  • returnDynMatrices::Bool: Flag to save linearized K/C/M matrices for the design

Outputs:

  • freq::Array{<:float}: sorted modal frequencies (Hz)
  • damp::Array{<:float}: sorted modal damping
  • imagCompSign::Array{<:float}: sign of imaginary component of eigenvalues
  • U_x_0::Array{<:float}: NnodesxNmodes in-phase mode shape x
  • U_y_0::Array{<:float}: NnodesxNmodes in-phase mode shape y
  • U_z_0::Array{<:float}: NnodesxNmodes in-phase mode shape z
  • theta_x_0::Array{<:float}: NnodesxNmodes in-phase mode shape rotation about x
  • theta_y_0::Array{<:float}: NnodesxNmodes in-phase mode shape rotation about y
  • theta_z_0::Array{<:float}: NnodesxNmodes in-phase mode shape rotation about z
  • U_x_90::Array{<:float}: NnodesxNmodes out-of-phase mode shape x
  • U_y_90::Array{<:float}: NnodesxNmodes out-of-phase mode shape y
  • U_z_90::Array{<:float}: NnodesxNmodes out-of-phase mode shape z
  • theta_x_90::Array{<:float}: NnodesxNmodes out-of-phase mode shape rotation about x
  • theta_y_90::Array{<:float}: NnodesxNmodes out-of-phase mode shape rotation about y
  • theta_z_90::Array{<:float}: NnodesxNmodes out-of-phase mode shape rotation about z
OWENSFEA.reducedOrderModelMethod

reducedOrderModel(elStorage,feamodel,mesh,el,displ)

This function executes a reduced order model analysis.

#Input

  • elStorage object containing stored element matrices
  • feamodel object containing feamodel information
  • mesh object containing mesh information
  • el object containing element information
  • displ displacement vector for use in pre-stressed analysis

#Output

  • rom object containing a reduced order feamodel
OWENSFEA.setInitialConditionsMethod
setInitialConditions(initCond,u,numDOFPerNode)

sets initial conditions

#Input

  • initCond: array containing initial conditions initCond(i,1) node number for init cond i initCond(i,2) local DOF number for init cond i initCond(i,3) value for init cond i
  • u: displacement vector for each dof
  • numDOFPerNode: number of degrees of freedom per node

#Output

  • u: displacement vector modified for initial conditions
OWENSFEA.setPrescribedConditionsMethod
prescribed_conditions = setPrescribedConditions(mesh;pBC=zeros(2,2),Fexternal=[],ForceDof=[])

Internal, maps OWENS boundary conditions and applied forces to the GXBeam PrescribedConditions input

#Input

  • mesh::FEAModel.mesh: Input turbine mesh
  • pBC::FEAModel.BC.pBC: Boundary conditions
  • Fexternal::Array{Float64}: Applied forces to the mesh, 1D array ordered node 1 Dof 1-6, node 2 Dof 1-2, etc.
  • ForceDof::Array{Float64}: Degrees of freedom aligned with Fexternal, currently is unused, assumes Fexternal uses the full DOF array.

#Output

  • prescribed_conditions::GXBeam.PrescribedConditions: the boundary conditions/applied forces used by GXBeam see ?GXBeam.PrescribedConditions
OWENSFEA.staticAnalysisMethod

staticAnalysis(feamodel,mesh,el,displ,Omega,OmegaStart,elStorage; reactionNodeNumber=1, OmegaDot=0.0, Fdof=[1], Fexternal=[0.0])

This function performs a static analysis and returns displacement values and a flag denoting successful/unsuccessful analysis

#Inputs

  • feamodel: object containing feamodel information
  • mesh: object containing mesh information
  • el: object containing element information
  • displ: displacement vector for use in pre-stressed analysis
  • Omega: rotor speed (Hz)
  • OmegaStart: rotor speed (Hz) from previous analysis if stepping through various rotor speeds, may be useful in load stepping
  • elStorage: previously calculated element system matrices
  • reactionNodeNumber::Int: optional, node at which to calculate reaction force
  • OmegaDot::Float: Steady State Rotational Acceleration
  • Fdof::Array{<:Int}: Global Dofs where Fexternal is acting, where max dof = nelem*ndof
  • Fexternal{<:Float}: Forces or moments associated with the Fdofs specified

#Outputs

  • displ: vector of displacemetns
  • staticAnalysisSuccessful: boolean flag denoting successful static analysis
OWENSFEA.structuralDynamicsTransientMethod
structuralDynamicsTransient(feamodel,mesh,el,dispData,Omega,OmegaDot,time,delta_t,elStorage,Fexternal,Fdof,CN2H,rbData)

performs unsteady structural dynamics analysis

Inputs

  • feamodel::: object containing feamodel data
  • mesh::: object containing mesh data
  • el::: object containing element data
  • dispData::: object containing displacement data
  • Omega::: rotor speed (Hz)
  • OmegaDot::: rotor acceleratin (Hz)
  • time::: current simulation time
  • delta_t::: time step size
  • elStorage::: object containing stored element data
  • Fexternal::: vector containing external force values
  • Fdof::: vector containing global DOF numbering associated with external force values
  • CN2H::: transformation matrix from inertial frame to hub frame
  • rbData::: vector containing rigid body displacement, velocity, and acceleration

Outputs:

  • elStrain::ElStrain: see ?ElStrain strain for element at end of time step
  • dispOut::DispOut: see ?DispOut displacement data at end of time step
  • FReaction_sp1::: vector containing reaction force at turbine base at end of time step
OWENSFEA.structuralDynamicsTransientROMMethod
structuralDynamicsTransientROM(feamodel,mesh,el,dispData,Omega,OmegaDot,time,delta_t,elStorage,rom,Fexternal,Fdof,CN2H,rbData)

Performs transient structural dynamics analysis using a reduced order feamodel (ROM).

#Input

  • feamodel: object containing feamodel data
  • mesh: object containing mesh data
  • el: object containing element data
  • dispData: object containing displacement data
  • Omega: rotor speed (Hz)
  • OmegaDot: rotor acceleratin (Hz)
  • time: current simulation time
  • delta_t: time step size
  • elStorage: object containing stored element data
  • rom: object containing reduced order feamodel represnetation
  • Fexternal: vector containing external force values
  • Fdof: vector containing global DOF numbering associated with external force values
  • CN2H: transformation matrix from inertial frame to hub frame
  • rbData: vector containing rigid body displacement, velocity, and acceleration

#Output

  • dispOut: object containing displacement data at end of time step
  • FReaction_sp1: vector containing reaction force at turbine base at end of time step
OWENSFEA.timeIntegrateSubSystemEffMethod

timeIntegrateSubSystemEff(M,K,C,F,timeInt,u,udot,uddot)

Performs integration of a system using the Newmark-Beta method(constant-average acceleration sceheme). The integration parameters are calculated before hand and store in the timeInt object.

#Input

  • M system mass matrix
  • K system sttiffness matrix
  • C system damping matrix
  • F system force vector
  • timeInt object containing time integraton parameters
  • u displacement at beginning of time step
  • udot velocity at beginning of time step
  • uddot acceleration at beginning of time step

#Output

  • unp1: displacement at end of time step
  • udotnp1: velocity at end of time step
  • uddotnp1: acceleration at end of time step
OWENSFEA.updateLoadStepMethod
updateLoadStep(iterationCount,loadStepParams,loadStep,loadStepPrev,loadStepCount,displCopy,displ)

Updates the load stepping parameter whether through means of adaptive loadstepping or a specified load step profile.

#Input

  • iterationCount number of iterations for current load step
  • loadStepParams struct containing load step parameters
  • loadStep load step value for current load step
  • loadStepPrev load step value for previous load st ep
  • loadStepCount number of load steps performed up to this point
  • displPrev converged displacement vector form previous load step
  • displ displacement vector at current load step

#Output

  • loadStep new load step value
  • loadStepPrev load step value for previous load step
  • displ most up to date displacement vector in load stepping procedure
  • displPrev displacement vector at end of previous load step
  • staticAnalysisSuccessful boolean flag, true if load step was completed successfully
  • staticAnalysisComplete boolean flag, true if analysis is complete