OWENSFEA
Types and functions
OWENSFEA.BC_struct
— TypeInternal, boundary condition data, see ?FEAModel for pBC
OWENSFEA.ConcNDL1D
— TypeInternal, NodalTerms node number, local dof (diagonal), and value
OWENSFEA.ConcNDL2D
— TypeInternal, NodalTerms node number, local dof1, local dof2, and value
OWENSFEA.DispData
— TypeInternal, displacement, velocity, and acceleration for each element
OWENSFEA.DispOut
— TypeDispOut(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.El
— TypeInternal, see ?Ort and ?SectionPropsArray
OWENSFEA.ElInput
— TypeInternal, inputs to Timoshenko element, inputs given by FEAModel struct and mesh
OWENSFEA.ElOutput
— TypeInternal, timoshenko element output matrices
OWENSFEA.ElStorage
— TypeInternal, Timoshenko element matrices
OWENSFEA.ElStrain
— TypeElStrain(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 directionepsilon_y::float
: epsilon_y strain in the y directionepsilon_z::float
: epsilon_z strain in the z directionkappa_x::float
: kappa_x curvature in the x directionkappa_y::float
: kappa_y curvature in the y directionkappa_z::float
: kappa_z curvature in the z direction
Outputs:
none
:
OWENSFEA.FEAModel
— MethodFEAModel(;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 valueaeroElasticOn::Bool
: Include simplified flutter calculataions in the timoshenko element matricesguessFreq::float
: aeroelastic starting guess, only used if aeroElasticOnairDensity::float
: working fluid densitygravityOn::Bool orArray{<:float}
: vector of 3 or flag to include distributed gravity acceleration (9.81m/s) in the negative z-directionnlOn::Bool
: flag for solver to calculate deflection induced stiffness changes and associated convergance to the coupled solutionspinUpOn::Bool
: flag to perform a static analysis (warm start) prior to performing modal analysisoutFilename::string
: /path/to/desired/output/filename if it doesn't exist already it is created, if exists, is overwrittenRayleighAlpha::float
: Rayleigh alpha damping used in timoshenko beam damping matrixRayleighBeta::float
: Rayleigh beta damping used in timoshenko beam damping matrixelementOrder::int
: order of element: 1 linear, 2 quadraticjoint::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 jointThetaplatformTurbineConnectionNodeNumber::int
: node at which reaction forces are calculatedjointTransform
: not used as an input, is calculated, local transform between dependent and active DOFs for nodes associated with a particular jointreducedDOFList::Array{<:int}
: not used as an input, is calculated, map of original DOF numbering to reduced DOF numberingnumDOFPerNode::int
: number of degrees of freedom per nodenumNodes::int
: total number of nodes in the meshnumModes::int
: number of modes to calculatenlParams::NlParams
: optional there in case the Nlparams struct is passed in, should be cleaned up since redundantalpha::float64
: optional newmark beta alpha parameter,If TD, use 0.25gamma::float64
: optional newmark beta gamma parameter, if static, use 0. If hydro, use 1.0pBC::Array{<:float}
: Nx3 array consisting of node, local dof, specified displacement value for the boundary conditionnodalTerms
: Concentrated nodal terms, should be replaced with the nodal input data array and the calc done internallyiterationType::string
: FEA displacement update calculation, Newton Raphson "NR", Direct Iteration "DI"adaptiveLoadSteppingFlag
: Unused, should be removedtolerance::float
: FEA total mesh unsteady analysis convergence tolerance for a timestepmaxIterations
: FEA total mesh unsteady analysis convergence max iterations for a timestepmaxNumLoadSteps
: used in static (steady state) analysisminLoadStepDelta
: used in static (steady state) analysisminLoadStep
: used in static (steady state) analysisprescribedLoadStep
: used in static (steady state) analysispredef::Bool
: will update the elStorage array if passed into Unsteady() with the nonlinear strain stiffening, to be used for subsequent analyses
Outputs:
none
:
OWENSFEA.Mesh
— TypeMesh(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 elementsnumNodes::int
: total number of nodesx::Array{<:float}
: Nodal x positiony::Array{<:float}
: Nodal y positionz::Array{<:float}
: Nodal z positionelNum::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-strutmeshSeg::Array{<:int}
: number of nodes within each segment, with segments consisting of tower, blade 1 2 etc, strutsstructuralSpanLocNorm::Array{<:float}
: Should be named heigh loc norm - unitized position along the blade height, used for aeroload mappingstructuralNodeNumbers::Array{<:int}
: Node numbers associated with blades for aero loads mappingstructuralElNumbers::Array{<:int}
: Element numbers associated with blades for aero loads mappingnonRotating::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 wireshubNodeNum::int
: Node number where the rotating part of the turbine starts, assumes meshing always starts with tower, then blades, etc.
Outputs:
none
:
OWENSFEA.NlParams
— TypeNlParams(iterationType,adaptiveLoadSteppingFlag,tolerance,maxIterations,maxNumLoadSteps,minLoadStepDelta,minLoadStep,prescribedLoadStep)
See ?FEAModel
OWENSFEA.NodalTerms
— TypeInternal, see ?FEAModel for NodalTerms
OWENSFEA.Ort
— TypeOrt(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 withOffset::Array{<:float}
: hub frame coordinate of node 1 of the element
Outputs:
none
:
OWENSFEA.ROM
— TypeInternal, ROM data
OWENSFEA.SectionPropsArray
— TypeSectionPropsArray(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 approximationtwist::Array{<:float}
: element twist (rad)rhoA::Array{<:float}
: rho * A in standard SI unitsEIyy::Array{<:float}
: E * IyyEIzz::Array{<:float}
: E * IzzGJ::Array{<:float}
: G * JEA::Array{<:float}
: E * ArhoIyy::Array{<:float}
: rho * IyyrhoIzz::Array{<:float}
: rho * IzzrhoJ::Array{<:float}
: rho * Jzcm::Array{<:float}
: z location of center of massycm::Array{<:float}
: y location of center of massa::Array{<:float}
: possibly lift slopeEIyz::Array{<:float}
: E * Iyzalpha1::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 codealpha2::Array{<:float}
: doesn't appear to be usedalpha3::Array{<:float}
: doesn't appear to be usedalpha4::Array{<:float}
: doesn't appear to be usedalpha5::Array{<:float}
: doesn't appear to be usedalpha6::Array{<:float}
: doesn't appear to be usedrhoIyz::Array{<:float}
: rho * Iyzb::Array{<:float}
: used in flutter approximation, possibly a chord or thickness valuea0::Array{<:float}
: zero lift angle of attack, used in flutter approximationaeroCenterOffset::Array{<:float}
: doesn't appear to be usedxaf::Array{<:float}
: x airfoil coordinates (to scale)yaf::Array{<:float}
: y airfoil coordinates (to scale)
Outputs:
none
:
OWENSFEA.TimeInt
— TypeInternal, time integration terms
OWENSFEA.ConcMassAssociatedWithElement
— MethodConcMassAssociatedWithElement(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 elementjoint
joint array for nodal termsnodalMassTerms
listing of concentrated nodal mass termsnodalStiffnessTerms
listing of concentrated nodal stiffness termsnodalLoads
listing of concentrated nodal loads terms
#Output
mass
array of concentrated mass associated with elementstiff
array of concentrated stiffness associated with elementload
array of concentrated loads associated with elementmodJoint
modified joint object removing nodal terms that have/will be applied to the element calculationsmodNodalMassTerms
modified nodal mass object removing nodal terms that have/will be applied to the element calculationsmodalStiffnessTerms
modified nodal stiffness object removing nodal terms that have/will be applied to the element calculationsmodNodalLoads
modified nodal loads object removing nodal terms that have/will be applied to the element calculations
OWENSFEA.ModalOutput
— MethodwriteOutput(freq,damp,phase1,phase2,imagComponentSign,fid)
Internal, writes an output file and or formats an output for modal analysis.
#Input
freq
: array of modal frequenciesdamp
: array of modal damping ratiosphase1
: array of in phase mode shapesphase2
: array of out of phase mode shapesimagComponentSign
: array of sign of imaginary componentsfid
: file identifier for output
#Output
freqSorted
: array of sorted(by frequency) modal frequenciesdampSorted
: array of sorted(by frequency) modal damping ratiosimagCompSignSorted
: array of sorted(by frequency) of imaginarycomponentSign arrayU_x_0
: see ?Modal outputsU_y_0
: see ?Modal outputsU_z_0
: see ?Modal outputstheta_x_0
: see ?Modal outputstheta_y_0
: see ?Modal outputstheta_z_0
: see ?Modal outputsU_x_90
: see ?Modal outputsU_y_90
: see ?Modal outputsU_z_90
: see ?Modal outputstheta_x_90
: see ?Modal outputstheta_y_90
: see ?Modal outputstheta_z_90
: see ?Modal outputs
OWENSFEA.adaptiveLoadStepping
— MethodInternal, performs updates a loadstep adaptively, see ?updateLoadStep
OWENSFEA.applyBC
— MethodapplyBC(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 matrixFg
assembled global load vectorBC
struct of boundary condition informationu
global displacement vectoriterationType
for nonlinear analysis, not used in BLASTnumDofPerNode
number of degrees of freedom per node
#Output
Kg
global stiffness matrix with boundary conditionsFg
global load vector with boundary condition
OWENSFEA.applyBCModal
— MethodapplyBCModal(K,BC,numDofPerNode)
Internal, applies boundary conditions to a system matrix for modal analysis
Inputs
K
: assembled global system matrixBC
: struct of boundary condition informationnumDofPerNode
: number of degrees of freedom per node
Outputs:
Knew
global system matrix with boundary conditions
OWENSFEA.applyConcentratedTerms
— MethodapplyConcentratedTerms(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 filenamedata
: 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.applyConstraints
— MethodInternal, this function transforms a matrix by the transformation matrix to enforce joint constraints
OWENSFEA.applyConstraintsVec
— MethodInternal, this function transforms a vector by the transformation matrix to enforce joint constraints
OWENSFEA.assembly!
— Methodassembly(Ke,Fe,conn,numNodesPerEl,numDOFPerNode,Kg,Fg)
Internal, assembles the element matrix and load vector into the global system of equations
#Input
Ke
: element matrixFe
: element vectorconn
: element connectivitynumNodesPerEl
: number of nodes per elementnumDofPerNode
: number of degrees of freedom per nodeKg
: global system matrixFg
: global load vector
#Output
Kg
: global system matrix with assembled elementFg
: global load vector with assembled element
OWENSFEA.assemblyMatrixOnly
— MethodassemblyMatrixOnly(Ke,conn,numNodesPerEl,numDOFPerNode,Kg)
Internal, assembles the element matrix into the global system of equations
Inputs
Ke
: element matrixconn
: element connectivitynumNodesPerEl
: number of nodes per elementnumDofPerNode
: number of degrees of freedom per nodeKg
: global system matrix
Outputs:
Kg
: global system matrix with assembled element
OWENSFEA.autoCampbellDiagram
— Methodfrequencies = 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 optionsmymesh::OWENSFEA.Mesh
: a previously generated turbine meshmyel::OWENSFEA.El
: the element properties associated with that meshsystem::GXBeam.System
: the converted GXBeam system from the mesh and elassembly::GXBeam.AssemblyState
: the converted GXBeam assembly from the mesh and elsections::Array{Float64, 3}
: the 3D point cloud to be converted to VTK formatminRPM::Float64
: minimum RPM to be run, e.x. 0.0maxRPM::Float64
: maximum RPM to be run e.x. 40.0NRPM::Int64
: number of linear discretizations of RPM e.x. 9 must be intvtksavename::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 intsaveRPM::Array{Int64}
: The RPMs to save in the VTK outputs e.x. [1,3,5] must be intmode_scaling::Float64
: The mode scaling in the VTK outputs e.x. 500.0
#Outputs
frequency::Array{Float64}
: The output modal frequencies
OWENSFEA.calcUnorm
— MethodThis function calculates a relative norm between two vectors: unew and uold
OWENSFEA.calculateBCMap
— MethodcalculateBCMap(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 conditionspBC
array of boundary condition datanumDofPerNode
number of degrees of freedom per nodereducedDofList
array of reduced DOF numbering
#Output
elStorage
map for boundary conditions between full and reduced dof list
OWENSFEA.calculateElement1!
— MethodInternal, general routine to calculate an element matrix
OWENSFEA.calculateElementMass
— MethodInternal, calculates element mass properties.
OWENSFEA.calculateLambda
— MethodcalculateLambda(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.calculateLoadVecFromDistForce
— MethodcalculateLoadVecFromDistForce(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-positionxloc::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 forceextDistF3Node::Array{<:float}
: turbine Normal forceextDistF4Node::Array{<:float}
: turbine M25 moment
#Output
Fe::Array{float}
: 6x1 Force on element in element FOR
OWENSFEA.calculateROM
— MethodcalculateROM(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 datamesh
object containing mesh datael
object containing elementdatadispl
displacement vectorrbData
: vector containing rigid body displacement, velocity, and accelerationelStorage
object containing stored element datacountedNodes
prevents applied nodal terms from double counting
#Output
rom
object containing reduced order model data
OWENSFEA.calculateROMGyric
— MethodcalculateROMGyric(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 datamesh
: object containing mesh datael
: object containing elementdatadispl
: displacement vectorrbData
: vector of hub frame accel (1-3), angular velocity components (4-6), and angular accleration (7-9)elStorage
: object containing stored element datarom0
: object containing parked/conventional reduced order feamodelcountedNodes
: prevents applied nodal terms from double counting
#Output
rom
: object containing reduced order feamodel data
OWENSFEA.calculateReactionForceAtNode
— MethodcalculateReactionForceAtNode(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 atmodel
: object containing model datamesh
: object containing mesh dataelStorage
: object containing stored element datael
: object containing element datatimeInt
: object containing time integration parametersdispData
: object containing displacement datadispl_iter
: converged displacement solutionrbData
: vector containing rigid body displacement, velocity, and accelerationOmega
: rotor speed (Hz)OmegaDot
: rotor acceleratin (Hz)CN2H
: transformation matrix from inertial frame to hub framecountedNodes
: prevents nodal terms from being double counted
#Output
cummulativeForce
: vector containing reaction force at nodeNum
OWENSFEA.calculateReducedDOFVector
— MethodInternal, searches over all DOFs in a structural model and determines and returns "dofVector" containing only unconstrained DOFs
OWENSFEA.calculateShapeFunctions
— MethodcalculateShapeFunctions(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 quadraticxi
guass point values to evaluate shape functions atx
nodal coordinates in physical length domain
#Output
N
shape function value at specified gauss pointsp_N_x
shape function derivative w.r.t physical length domain at specified gauss pointsJac
Jacobian for mat between local element domain and physical length domain.
OWENSFEA.calculateStrainForElements
— MethodInternal calculates element strains
OWENSFEA.calculateStructureMassProps
— MethodcalculateStructureMassProps(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 structurestructureMOI::float
moment of inertia tensor of structgurestructureMassCenter::float
center of mass of structure
OWENSFEA.calculateTimoshenkoElementInitialRun
— MethodcalculateTimoshenkoElementInitialRun(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.calculateTimoshenkoElementNL
— MethodcalculateTimoshenkoElementNL(input,elStorage;predef=nothing)
Internal, performs nonlinear element calculations.
#Inputs
input::ElInput
: see ?ElInputelStorage::ElStorage
: see ?ElStoragepredef::Bool
: optional, if true, mutates ElStorage to include the nonlinear strain stiffening
#Outputs
ElOutput
: see ?ElOutput
OWENSFEA.calculateTimoshenkoElementNLSS
— MethodcalculateTimoshenkoElementNLSS(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.calculateTimoshenkoElementStrain
— MethodcalculateTimoshenkoElementStrain(elementOrder,nlOn,xloc,sectionProps,sweepAngle,coneAngle,rollAngle,aeroSweepAngle,disp)
Internal, calculates element strain for a Timoshenko element
#Outputs
ElStrain
: See ?ElStrain
OWENSFEA.calculateVec1!
— MethodInternal, general routine to calculate an element vector
OWENSFEA.constructReducedDispVecFromEigVec
— MethodconstructReducedDispVecFromEigVec(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.constructReducedDispVectorMap
— MethodInternal, creates a map of unconstrained DOFs between a full listing and reduced listing (after constraints have been applied)
OWENSFEA.createJointTransform
— MethodcreateJointTransform(joint,numNodes,numDofPerNode)
Internal, calculates the JointTransform of a structural system.
#Input
joint
: object containing joint datanumModes
: number of nodes in meshnumDofPerNode
: number of degrees of freedom per node
#Output
jointTransform
: joint transformation matrixreducedDOF
: map of original DOF numbering to reduced DOF numbering
OWENSFEA.createTda
— MethodInternal, 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.determineActiveDofsFromSlaveNode
— MethodInternal, determines the local master DOF associated with a local slave DOF.
OWENSFEA.extractFreqDamp
— MethodextractFreqDamp(val,vec,numDOFPerNode,jointTransform,reducedDOFList,BC,analysisType)
Internal, calculates the eigenvalues and vectors of a structural dynamic system
Inputs
val
: eigenvaluevec
: eigenvectornumDOFPerNode
: number of degrees of freedom per nodejointTransform
: joint transformation matrix from reduced to full DOF listreducedDOFList
: listing of reduced DOFsBC
: boundary condition object containing boundary condition infoanalysisType
: analysis type
Outputs:
freq
: modal frequencydamp
: modal dampingphase1
: 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.extractdaInfo
— MethodInternal, 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.findElementsAssociatedWithNodeNumber
— MethodfindElementsAssociatedWithNodeNumber(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 atconn
object containing mesh connectivityjointData
object containing joint information
#Output
elList
array containing a list of element numbers associated with nodeNumlocalNode
array containing the local node number that correspond to nodeNum in the list of associated elements
OWENSFEA.getElementConcTerms!
— MethodInternal, gets the concentrated terms without double counting
OWENSFEA.getGP
— MethodgetGP(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 frameweight
: associated weights for quad point coordinate
OWENSFEA.getNodeMaps
— MethodInternal, gets node mapping
OWENSFEA.initialElementCalculations
— MethodinitialElementCalculations(feamodel,el,mesh)
performs initial element calculation for use later in analysis for efficiency gains.
Inputs
feamodel::FEAmodel
: see ?Feamodelel::El
: see ?Elmesh::Mesh
: see ?Mesh
Outputs:
elStorage:ElStorage
: see ?ElStorage
OWENSFEA.interpolateVal
— MethodInternal, linear interpolation
OWENSFEA.linearAnalysisModal
— MethodInternal, see ?modal
OWENSFEA.makeBCdata
— MethodmakeBCdata(pBC,numNodes,numDofPerNode,reducedDOFList,jointTransform)
Internal, usese the pBC matrix and calculates/stores boundary condition data
#Input
pBC
See ?FEAModel.pBCnumNodes
number of nodes in structural modelnumDofPerNode
number of degrees of freedom per nodereducedDOFList
joint transformation matrix from reduced to full DOF listjointTransform
listing of reduced DOFs
#Output
BC:BC_struct
see ?BC_struct
OWENSFEA.mapMatrixNonSym
— MethodInternal, function to form total stifness matrix and transform to desired DOF mapping
OWENSFEA.mapMatrixNonSym2
— MethodInternal, function to form total stifness matrix and transform to desired DOF mapping
OWENSFEA.mapVector
— MethodInternal, forms total force vector and transform to desired DOF mapping
OWENSFEA.modal
— Methodmodal(feamodel,mesh,el;Omega=0.0,displ=zeros(mesh.numNodes*6),OmegaStart=0.0,returnDynMatrices=false)
Modal analysis
Inputs
feamodel::FEAModel
: see ?FEAModelmesh::Mesh
: see ?Meshel::El
: see ?ElOmega::float
: Rotational rate in Hzdispl::Array{<:float}
: zeros(mesh.numNodes*6) initial (warm start) displacements for each dofOmegaStart::float
: rotor speed (Hz) from previous analysis if stepping through various rotor speeds, may be useful in load steppingreturnDynMatrices::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 dampingimagCompSign::Array{<:float}
: sign of imaginary component of eigenvaluesU_x_0::Array{<:float}
: NnodesxNmodes in-phase mode shape xU_y_0::Array{<:float}
: NnodesxNmodes in-phase mode shape yU_z_0::Array{<:float}
: NnodesxNmodes in-phase mode shape ztheta_x_0::Array{<:float}
: NnodesxNmodes in-phase mode shape rotation about xtheta_y_0::Array{<:float}
: NnodesxNmodes in-phase mode shape rotation about ytheta_z_0::Array{<:float}
: NnodesxNmodes in-phase mode shape rotation about zU_x_90::Array{<:float}
: NnodesxNmodes out-of-phase mode shape xU_y_90::Array{<:float}
: NnodesxNmodes out-of-phase mode shape yU_z_90::Array{<:float}
: NnodesxNmodes out-of-phase mode shape ztheta_x_90::Array{<:float}
: NnodesxNmodes out-of-phase mode shape rotation about xtheta_y_90::Array{<:float}
: NnodesxNmodes out-of-phase mode shape rotation about ytheta_z_90::Array{<:float}
: NnodesxNmodes out-of-phase mode shape rotation about z
OWENSFEA.reducedOrderModel
— MethodreducedOrderModel(elStorage,feamodel,mesh,el,displ)
This function executes a reduced order model analysis.
#Input
elStorage
object containing stored element matricesfeamodel
object containing feamodel informationmesh
object containing mesh informationel
object containing element informationdispl
displacement vector for use in pre-stressed analysis
#Output
rom
object containing a reduced order feamodel
OWENSFEA.setInitialConditions
— MethodsetInitialConditions(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 iu
: displacement vector for each dofnumDOFPerNode
: number of degrees of freedom per node
#Output
u
: displacement vector modified for initial conditions
OWENSFEA.setPrescribedConditions
— Methodprescribed_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 meshpBC::FEAModel.BC.pBC
: Boundary conditionsFexternal::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.staticAnalysis
— MethodstaticAnalysis(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 informationmesh
: object containing mesh informationel
: object containing element informationdispl
: displacement vector for use in pre-stressed analysisOmega
: rotor speed (Hz)OmegaStart
: rotor speed (Hz) from previous analysis if stepping through various rotor speeds, may be useful in load steppingelStorage
: previously calculated element system matricesreactionNodeNumber::Int
: optional, node at which to calculate reaction forceOmegaDot::Float
: Steady State Rotational AccelerationFdof::Array{<:Int}
: Global Dofs where Fexternal is acting, where max dof = nelem*ndofFexternal{<:Float}
: Forces or moments associated with the Fdofs specified
#Outputs
displ
: vector of displacemetnsstaticAnalysisSuccessful
: boolean flag denoting successful static analysis
OWENSFEA.structuralDynamicsTransient
— MethodstructuralDynamicsTransient(feamodel,mesh,el,dispData,Omega,OmegaDot,time,delta_t,elStorage,Fexternal,Fdof,CN2H,rbData)
performs unsteady structural dynamics analysis
Inputs
feamodel::
: object containing feamodel datamesh::
: object containing mesh datael::
: object containing element datadispData::
: object containing displacement dataOmega::
: rotor speed (Hz)OmegaDot::
: rotor acceleratin (Hz)time::
: current simulation timedelta_t::
: time step sizeelStorage::
: object containing stored element dataFexternal::
: vector containing external force valuesFdof::
: vector containing global DOF numbering associated with external force valuesCN2H::
: transformation matrix from inertial frame to hub framerbData::
: vector containing rigid body displacement, velocity, and acceleration
Outputs:
elStrain::ElStrain
: see ?ElStrain strain for element at end of time stepdispOut::DispOut
: see ?DispOut displacement data at end of time stepFReaction_sp1::
: vector containing reaction force at turbine base at end of time step
OWENSFEA.structuralDynamicsTransientROM
— MethodstructuralDynamicsTransientROM(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 datamesh
: object containing mesh datael
: object containing element datadispData
: object containing displacement dataOmega
: rotor speed (Hz)OmegaDot
: rotor acceleratin (Hz)time
: current simulation timedelta_t
: time step sizeelStorage
: object containing stored element datarom
: object containing reduced order feamodel represnetationFexternal
: vector containing external force valuesFdof
: vector containing global DOF numbering associated with external force valuesCN2H
: transformation matrix from inertial frame to hub framerbData
: vector containing rigid body displacement, velocity, and acceleration
#Output
dispOut
: object containing displacement data at end of time stepFReaction_sp1
: vector containing reaction force at turbine base at end of time step
OWENSFEA.timeIntegrateSubSystemEff
— MethodtimeIntegrateSubSystemEff(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 matrixK
system sttiffness matrixC
system damping matrixF
system force vectortimeInt
object containing time integraton parametersu
displacement at beginning of time stepudot
velocity at beginning of time stepuddot
acceleration at beginning of time step
#Output
unp1
: displacement at end of time stepudotnp1
: velocity at end of time stepuddotnp1
: acceleration at end of time step
OWENSFEA.updateLoadStep
— MethodupdateLoadStep(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 steploadStepParams
struct containing load step parametersloadStep
load step value for current load steploadStepPrev
load step value for previous load st eploadStepCount
number of load steps performed up to this pointdisplPrev
converged displacement vector form previous load stepdispl
displacement vector at current load step
#Output
loadStep
new load step valueloadStepPrev
load step value for previous load stepdispl
most up to date displacement vector in load stepping proceduredisplPrev
displacement vector at end of previous load stepstaticAnalysisSuccessful
boolean flag, true if load step was completed successfullystaticAnalysisComplete
boolean flag, true if analysis is complete