Detailed Inputs

In this example, we show the second level of what is going on behind the precompiled binary This would be appropriate if you need more customization in the run and design parameters than the input file currently allows, but your design still fits within the setupOWENS helper function etc.

Tip

This example is also available as a Jupyter notebook todo: get link working:

First we import the packages. If "using" was employed, then all of the functions of the packages specified would be made available, but "import" requires PackageName.FunctionName to be used unless the function was explicitely exported in the package. If "include("filepath/filename.jl)" is used, that is the same as copying and pasting. Please see the respective page on YAML input (TODO) for a description of the YAML inputs

import OWENS
import OWENSAero
#### import PyPlot
runpath = path = "/home/runner/work/OWENS.jl/OWENS.jl/examples/literate" # to run locally, change to splitdir(@__FILE__)[1]
# runpath = path = splitdir(@__FILE__)[1]
Inp = OWENS.MasterInput("$runpath/sampleOWENS.yml")

nothing

Unpack inputs, or you could directly input them here and bypass the file

verbosity = 1

analysisType = Inp.analysisType
turbineType = Inp.turbineType
eta = Inp.eta
Nbld = Inp.Nbld
towerHeight = Inp.towerHeight
rho = Inp.rho
Vinf = Inp.Vinf
controlStrategy = Inp.controlStrategy
RPM = Inp.RPM
Nslices = Inp.Nslices
ntheta = Inp.ntheta
structuralModel = Inp.structuralModel
ntelem = Inp.ntelem
nbelem = Inp.nbelem
ncelem = Inp.ncelem
nselem = Inp.nselem
ifw = Inp.ifw
WindType = Inp.WindType
AModel = Inp.AModel
windINPfilename = "$(path)$(Inp.windINPfilename)"
ifw_libfile = Inp.ifw_libfile
if ifw_libfile == "nothing"
    ifw_libfile = nothing
end
Blade_Height = Inp.Blade_Height
Blade_Radius = Inp.Blade_Radius
numTS = Inp.numTS
delta_t = Inp.delta_t
NuMad_geom_xlscsv_file_twr = "$(path)$(Inp.NuMad_geom_xlscsv_file_twr)"
NuMad_mat_xlscsv_file_twr = "$(path)$(Inp.NuMad_mat_xlscsv_file_twr)"
NuMad_geom_xlscsv_file_bld = "$(path)$(Inp.NuMad_geom_xlscsv_file_bld)"
NuMad_mat_xlscsv_file_bld = "$(path)$(Inp.NuMad_mat_xlscsv_file_bld)"
NuMad_geom_xlscsv_file_strut = "$(path)$(Inp.NuMad_geom_xlscsv_file_strut)"
NuMad_mat_xlscsv_file_strut = "$(path)$(Inp.NuMad_mat_xlscsv_file_strut)"
adi_lib = Inp.adi_lib
if adi_lib == "nothing"
    adi_lib = nothing
end
adi_rootname = "$(path)$(Inp.adi_rootname)"

B = Nbld
R = Blade_Radius#177.2022*0.3048 #m
H = Blade_Height#1.02*R*2 #m

shapeZ = collect(LinRange(0,H,Nslices+1))
shapeX = R.*(1.0.-4.0.*(shapeZ/H.-.5).^2)#shapeX_spline(shapeZ)

nothing

Call the helper function that builds the mesh, calculates the sectional properties, and aligns the sectional properties to the mesh elements,

mymesh,myel,myort,myjoint,sectionPropsArray,mass_twr, mass_bld,
stiff_twr, stiff_bld,bld_precompinput,
bld_precompoutput,plyprops_bld,numadIn_bld,lam_U_bld,lam_L_bld,
twr_precompinput,twr_precompoutput,plyprops_twr,numadIn_twr,lam_U_twr,lam_L_twr,aeroForces,deformAero,
mass_breakout_blds,mass_breakout_twr,system,assembly,sections,AD15bldNdIdxRng, AD15bldElIdxRng = OWENS.setupOWENS(OWENSAero,path;
    rho,
    Nslices,
    ntheta,
    RPM,
    Vinf,
    eta,
    B,
    H,
    R,
    shapeZ,
    shapeX,
    ifw,
    WindType,
    delta_t,
    numTS,
    adi_lib,
    adi_rootname,
    windINPfilename,
    ifw_libfile,
    NuMad_geom_xlscsv_file_twr,# = "$path/data/NuMAD_Geom_SNL_5MW_ARCUS_Cables.csv",
    NuMad_mat_xlscsv_file_twr,# = "$path/data/NuMAD_Materials_SNL_5MW_D_TaperedTower.csv",
    NuMad_geom_xlscsv_file_bld,# = "$path/data/NuMAD_Geom_SNL_5MW_ARCUS.csv",
    NuMad_mat_xlscsv_file_bld,# = "$path/data/NuMAD_Materials_SNL_5MW_D_Carbon_LCDT_ThickFoils_ThinSkin.csv",
    NuMad_geom_xlscsv_file_strut,
    NuMad_mat_xlscsv_file_strut,
    Htwr_base=towerHeight,
    strut_twr_mountpoint = [0.2,0.8],
    strut_bld_mountpoint = [0.2,0.8],
    ntelem,
    nbelem,
    ncelem,
    nselem,
    joint_type = 0,
    c_mount_ratio = 0.05,
    AModel, #AD, DMS, AC
    DSModel="BV",
    RPI=true,
    cables_connected_to_blade_base = true,
    meshtype = turbineType)

nothing
┌ Warning: Mesh warning: two points are directly on top of one another, consider adjusting number of elements to space out the mesh
└ @ OWENS ~/work/OWENS.jl/OWENS.jl/src/meshing_utilities.jl:392
┌ Warning: Mesh warning: two points are directly on top of one another, consider adjusting number of elements to space out the mesh
└ @ OWENS ~/work/OWENS.jl/OWENS.jl/src/meshing_utilities.jl:392
┌ Warning: Data for SN curve control points not found in material file columns 23:28 for stress in Mpa, 29:33 for cycles in log10
└ @ OWENS ~/work/OWENS.jl/OWENS.jl/src/fileio.jl:648
┌ Warning: Data for SN curve control points not found in material file columns 23:28 for stress in Mpa, 29:33 for cycles in log10
└ @ OWENS ~/work/OWENS.jl/OWENS.jl/src/fileio.jl:648
┌ Warning: Data for SN curve control points not found in material file columns 23:28 for stress in Mpa, 29:33 for cycles in log10
└ @ OWENS ~/work/OWENS.jl/OWENS.jl/src/fileio.jl:648
┌ Warning: Data for SN curve control points not found in material file columns 23:28 for stress in Mpa, 29:33 for cycles in log10
└ @ OWENS ~/work/OWENS.jl/OWENS.jl/src/fileio.jl:648
Opening AeroDyn-Inflow library at: /home/runner/.julia/packages/OWENSOpenFASTWrappers/OOzCv/src/../deps/openfast/build/modules/aerodyn/libaerodyn_inflow_c_binding

 **************************************************************************************************
 AeroDyn-Inflow library

 Copyright (C) 2024 National Renewable Energy Laboratory
 Copyright (C) 2024 Envision Energy USA LTD

 This program is licensed under Apache License Version 2.0 and comes with ABSOLUTELY NO WARRANTY.
 See the "LICENSE" file distributed with this software for details.
 **************************************************************************************************

 AeroDyn-Inflow library--128-NOTFOUND
 Compile Info:
  - Compiler: GCC version 11.4.0
  - Architecture: 64 bit
  - Precision: double
  - OpenMP: Yes, number of threads: 4/4
  - Date: Oct 17 2024
  - Time: 21:57:44
 Execution Info:
  - Date: 10/17/2024
  - Time: 22:14:31+0000

 Running ADI.
 Running AeroDyn.
 Running OLAF.
  - Directory:         /home/runner/work/OWENS.jl/OWENS.jl/docs/build/examples
  - RootName:          /home/runner/work/OWENS.jl/OWENS.jl/examples/literate/ExampleB.
  - Reading advanced options for OLAF:
 [WARN] Line ignored: "VERT2"   DEFAULT      -10       30.   400     5.      5.     1     5.136    15.136  100
  - OLAF regularization parameters (for wing 1):
    WingReg (min/max) :  12.2425 25.2472
    WakeReg (min/max) :  12.2425 25.2472
    k = alpha delta nu:   0.0368
 Running InflowWind.

Optionally, we can run the finite element solver with gemetrically exact beam theory via GXBeam.jl this requires that the OWENS style inputs are converted to the GXBeam inputs. This interface also includes the ability to output VTK files, which can be viewed in paraview. We have adapted this interface to work with OWENS inputs as well.

nothing

If the sectional properties material files includes cost information, that is combined with the density to estimate the overall material cost of of materials in the blades

if verbosity>0

    println("\nBlades' Mass Breakout")
    for (i,name) in enumerate(plyprops_bld.names)
        println("$name $(mass_breakout_blds[i]) kg, $(plyprops_bld.costs[i]) \$/kg: \$$(mass_breakout_blds[i]*plyprops_bld.costs[i])")
    end

    println("\nTower Mass Breakout")
    for (i,name) in enumerate(plyprops_twr.names)
        println("$name $(mass_breakout_twr[i]) kg, $(plyprops_twr.costs[i]) \$/kg: \$$(mass_breakout_twr[i]*plyprops_twr.costs[i])")
    end

    println("Total Material Cost Blades: \$$(sum(mass_breakout_blds.*plyprops_bld.costs))")
    println("Total Material Cost Tower: \$$(sum(mass_breakout_twr.*plyprops_twr.costs))")
    println("Total Material Cost: \$$(sum(mass_breakout_blds.*plyprops_bld.costs)+ sum(mass_breakout_twr.*plyprops_twr.costs))")

end

nothing

Blades' Mass Breakout
CLA_5500 20670.29648289118 kg, 2.06 $/kg: $42580.81075475583
CBX_2400 9730.417758654583 kg, 2.1 $/kg: $20433.877293174624
ETLX_2400 0.0 kg, 2.21 $/kg: $0.0
Airex_C70_55 789.5974763487575 kg, 7.23 $/kg: $5708.7897540015165
EBX_2400_x10 0.0 kg, 2.06 $/kg: $0.0
ETLX_2400_x10 0.0 kg, 2.1 $/kg: $0.0
Airex_C70_55_x10 0.0 kg, 7.23 $/kg: $0.0

Tower Mass Breakout
CLA_5500 4992.974575844168 kg, 2.06 $/kg: $10285.527626238985
CBX_2400 1217.5626883060443 kg, 2.1 $/kg: $2556.881645442693
ETLX_2400 0.0 kg, 2.21 $/kg: $0.0
Airex_C70_55 0.0 kg, 7.23 $/kg: $0.0
EBX_2400_x10 0.0 kg, 2.06 $/kg: $0.0
ETLX_2400_x10 0.0 kg, 2.1 $/kg: $0.0
Airex_C70_55_x10 0.0 kg, 7.23 $/kg: $0.0
Total Material Cost Blades: $68723.47780193196
Total Material Cost Tower: $12842.409271681678
Total Material Cost: $81565.88707361365

Here we apply the boundary conditions. For this case, with a regular cantelever tower, the tower base node which is 1 is constrained in all 6 degrees of freedom to have a displacement of 0. You can change this displacement to allow for things like pretension, and you can apply boundary conditions to any node.

pBC = [1 1 0
1 2 0
1 3 0
1 4 0
1 5 0
1 6 0]

nothing

There are inputs for the overall coupled simulation, please see the api reference for specifics on all the options

if AModel=="AD"
    AD15On = true
else
    AD15On = false
end

inputs = OWENS.Inputs(;analysisType = structuralModel,
tocp = [0.0,100000.1],
Omegaocp = [RPM,RPM] ./ 60,
tocp_Vinf = [0.0,100000.1],
Vinfocp = [Vinf,Vinf],
numTS,
delta_t,
AD15On,
aeroLoadsOn = 2)

nothing

Then there are inputs for the finite element models, also, please see the api reference for specifics on the options (TODO: ensure that this is propogated to the docs)

feamodel = OWENS.FEAModel(;analysisType = structuralModel,
outFilename = "none",
joint = myjoint,
platformTurbineConnectionNodeNumber = 1,
pBC,
nlOn = true,
numNodes = mymesh.numNodes,
RayleighAlpha = 0.05,
RayleighBeta = 0.05,
iterationType = "DI")

nothing

Here is where we actually call the unsteady simulation and where owens pulls the aero and structural solutions together and propogates things in time.

println("Running Unsteady")
t, aziHist,OmegaHist,OmegaDotHist,gbHist,gbDotHist,gbDotDotHist,FReactionHist,
FTwrBsHist,genTorque,genPower,torqueDriveShaft,uHist,uHist_prp,epsilon_x_hist,epsilon_y_hist,
epsilon_z_hist,kappa_x_hist,kappa_y_hist,kappa_z_hist = OWENS.Unsteady_Land(inputs;system,assembly,
topModel=feamodel,topMesh=mymesh,topEl=myel,aero=aeroForces,deformAero)

if AD15On #TODO: move this into the run functions
    OWENS.OWENSOpenFASTWrappers.endTurb()
end

nothing
Running Unsteady
Running in specified rotor speed mode

Simulation Time: 0.0 seconds of 0.49 seconds
 [INFO] FVW: Update States: reevaluation at the same starting time.  This will not print on
 subsequent occurences.

Simulation Time: 0.1 seconds of 0.49 seconds

Simulation Time: 0.2 seconds of 0.49 seconds

Simulation Time: 0.3 seconds of 0.49 seconds

Simulation Time: 0.4 seconds of 0.49 seconds
Simulation Complete.
 >>> FINAL WRITE

Like described above, we can output vtk files viewable in paraview. Here it is done for each time step and shows the deformations. Additionaly, there is a method to input custom values and have them show up on the vtk surface mesh for example, strain, or reaction force, etc. This is described in more detail in the api reference for the function and: TODO

println("Saving VTK time domain files")
OWENS.OWENSFEA_VTK("$path/vtk/SNLARCUS5MW_timedomain_TNBnltrue",t,uHist,system,assembly,sections;scaling=1,azi=aziHist)

nothing
Saving VTK time domain files

This helper function looks through all the loads and picks out the worst case safety factor in each of the stacks of composite lamina it also calculates analytical simply supported buckling safety factors

##########################################
#### Ultimate Failure #####
##########################################

massOwens,stress_U,SF_ult_U,SF_buck_U,stress_L,SF_ult_L,SF_buck_L,stress_TU,SF_ult_TU,
SF_buck_TU,stress_TL,SF_ult_TL,SF_buck_TL,topstrainout_blade_U,topstrainout_blade_L,
topstrainout_tower_U,topstrainout_tower_LtopDamage_blade_U,
topDamage_blade_L,topDamage_tower_U,topDamage_tower_L = OWENS.extractSF(bld_precompinput,
bld_precompoutput,plyprops_bld,numadIn_bld,lam_U_bld,lam_L_bld,
twr_precompinput,twr_precompoutput,plyprops_twr,numadIn_twr,lam_U_twr,lam_L_twr,
mymesh,myel,myort,Nbld,epsilon_x_hist,kappa_y_hist,kappa_z_hist,epsilon_z_hist,
kappa_x_hist,epsilon_y_hist;verbosity, #Verbosity 0:no printing, 1: summary, 2: summary and spanwise worst safety factor # epsilon_x_hist_1,kappa_y_hist_1,kappa_z_hist_1,epsilon_z_hist_1,kappa_x_hist_1,epsilon_y_hist_1,
LE_U_idx=1,TE_U_idx=6,SparCapU_idx=3,ForePanelU_idx=2,AftPanelU_idx=5,
LE_L_idx=1,TE_L_idx=6,SparCapL_idx=3,ForePanelL_idx=2,AftPanelL_idx=5,
Twr_LE_U_idx=1,Twr_LE_L_idx=1,
AD15bldNdIdxRng,AD15bldElIdxRng,strut_precompoutput=nothing) #TODO: add in ability to have material safety factors and load safety factors

nothing
Composite Ultimate and Buckling Safety Factors


UPPER BLADE SURFACE

Minimum Safety Factor on Surface: 0.30004180102744155
At time 0.034s at composite station 5 of 21 at lam 6 of 6
Maximum Damage per hr: 0.5142857152446287
At composite station 14 of 21 at lam 6 of 6


LOWER BLADE SURFACE

Minimum Safety Factor on Surface: 0.2998487967731362
At time 0.034s at composite station 5 of 21 at lam 6 of 6
Maximum Damage per hr: 0.7346938775718813
At composite station 9 of 21 at lam 4 of 6


UPPER TOWER

Minimum Safety Factor on tower Surface: 4.466109773579
At time 0.049s at composite station 7 of 21 at lam 1 of 1
Maximum Damage per hr: 0.22040816341419076
At composite station 21 of 21 at lam 1 of 1


Lower TOWER

Minimum Safety Factor on tower Surface: 4.663274990042024
At time 0.04s at composite station 1 of 21 at lam 1 of 1
Maximum Damage per hr: 0.1469387755538469
At composite station 13 of 21 at lam 1 of 1

This page was generated using Literate.jl.