DuctAPE.post_process
— Functionpost_process(
solver_options,
converged_states,
prepost_containers,
solve_container_caching,
solve_parameter_cache_vector,
solve_parameter_cache_dims,
operating_point,
reference_parameters,
boundary_layer_options,
A_bb_LU,
airfoils,
idmaps,
problem_dimensions,
multipoint_index;
write_outputs=options.write_outputs,
outfile=options.outfile,
checkoutfileexists=options.checkoutfileexists,
output_tuple_name=options.output_tuple_name,
verbose=options.verbose,
)
Post-process a converged nonlinear solve solution.
Arguments
solver_options::SolverOptionsType
: A SolverOptionsType object (also used for dispatch)converged_states::Vector{Float}
: the converged state variablesprepost_containers::NamedTuple
: the named tuple containing pre-allocated containers for the pre- and post-processing intermediate calculationssolve_container_cache::NamedTuple
: the cache and dimensions for intermediate values in the residual calculationsolve_parameter_cache_vector::Vector{Float}
: the applicably typed cache vector for the solve parameterssolve_parameter_cache_dims::NamedTuple
: the dimensions of the solver parametersoperating_point::OperatingPoint
: the operating point being analyzedreference_parameters::ReferenceParameters
: a ReferenceParameters objectBoundaryLayerOptions::BoundaryLayerOptions
: a BoundaryLayerOptions objectA_bb_LU::LinearAlgebra.LU
: LinearAlgebra LU factorization of the LHS matrixairfoils::Vector{AFType}
: A matrix of airfoil types associated with each of the blade elementsidmaps::NamedTuple
: A named tuple containing index mapping used in bookkeeping throughout solve and post-processproblem_dimensions::ProblemDimensions
: A ProblemDimensions object
Keyword Arguments
multipoint_index::Vector{Int}
: a one-dimensional vector containing the index of which multipoint analysis operating point is being analyzed.write_outputs=options.write_outputs::Vector{Bool}
: a vector with the same length as number of multipoints indicating if the outputs should be saved.outfile=options.outfile::Vector{String}
: a vector of file paths/names for where outputs should be writtencheckoutfileexists=options.checkoutfileexists::Bool
: a flag for whether existing files should be checked for or if blind overwriting is okay.output_tuple_name=options.output_tuple_name::Vector{String}
: the variable name(s) of the named tuple of outputs to be written.verbose::Bool=false
: flag to print verbose statements
Returns
outs::NamedTuple
: A named tuple containing all the output values including
bodies
panel_strengths
inviscid_thrust
viscous_drag
thrust_comp
total_thrust
induced_efficiency
cp_in
cp_out
cp_casing_in
cp_casing_out
casing_zpts
cp_nacelle_in
cp_nacelle_out
nacelle_zpts
cp_centerbody_in
cp_centerbody_out
centerbody_zpts
Vtot_in
Vtot_out
Vtot_prejump
vtot_body
vtot_jump
vtot_wake
vtot_rotors
Vtan_in
Vtan_out
vtan_casing_in
vtan_casing_out
vtan_nacelle_in
vtan_nacelle_out
vtan_centerbody_in
vtan_centerbody_out
boundary_layers
stagnation_indices
upper_initial_states
upper_solved_states
upper_solved_steps
lower_initial_states
lower_solved_states
lower_solved_steps
surface_length_upper
surface_length_lower
stag_point
split_ratio
separation_point_ratio_upper
separation_point_ratio_lower
cdc_upper
cdc_upper`cdc_lower
cdc_lower`vtdotpv
boundary_layer_functions_lower
boundary_layer_functions_upper
rotors
circulation
panel_strengths
efficiency
inviscid_thrust
inviscid_thrust_dist
viscous_thrust
viscous_thrust_dist
thrust
CT
inviscid_torque
inviscid_torque_dist
viscous_torque
viscous_torque_dist
torque
CQ
inviscid_power
inviscid_power_dist
viscous_power
viscous_power_dist
power
CP
cl
cd
alpha
beta1
blade_normal_force_per_unit_span
blade_tangential_force_per_unit_span
wake
panel_strengths
totals
thrust
torque
power
CT
CQ
CP
total_efficiency
ideal_efficiency
intermediate_solve_values
vz_rotor
vtheta_rotor
Cm_wake
reynolds
mach
Cz_rotor
Ctheta_rotor
Cmag_rotor
Gamma_tilde
H_tilde
deltaGamma2
deltaH
vz_wake
vr_wake
Cm_avg
reference_values
Vinf
Vref
Velocities
DuctAPE.get_body_tangential_velocities
— Functionget_body_tangential_velocities(
gamb,
gamw,
sigr,
ivb,
Vinf,
totnode,
totpanel,
nnode,
npanel,
tangent,
controlpoints,
endpanelidxs,
wake_panel_ids_along_centerbody_wake_interface,
wake_panel_ids_along_casing_wake_interface,
centerbody_panel_ids_along_centerbody_wake_interface,
duct_panel_ids_along_casing_wake_interface,
num_casing_panels,
)
Get the tangential velocities along the body surfaces.
Arguments
gamb::Vector{Float}
: the body panel strengthsgamw::Vector{Float}
: the wake panel strengthssigr::Vector{Float}
: the rotor panel strengthsivb::NamedTuple
: the unit induced velocities on the bodiesVinf::Vector{Float}
: one element vector containing the freestream magnitudetotnode::Int
: total number of nodes between all bodiestotpanel::Int
: total number of panels between all bodiesnnode::Vector{Int}
: number of nodes in each bodynpanel::Vector{Int}
: number of panels in each body.tangent::Matrix{Float}
: unit tangent vectors for each panelcontrolpoints::Matrix{Float}
: control point locations for each panelendpanelidxs::Matrix{Int}
: the indices of the first and last panels for each bodywake_panel_ids_along_centerbody_wake_interface::Vector{Int}
: the indices of the wake panels coincident with the centerbody panelswake_panel_ids_along_casing_wake_interface::Vector{Int}
: the indices of the wake panels coincident with the duct casing (inner surface) panelscenterbody_panel_ids_along_centerbody_wake_interface::Vector{Int}
: the indices of the centerbody panels coincident with the wake panelsduct_panel_ids_along_casing_wake_interface::Vector{Int}
: the indices of the duct panels coincident with the wake panelsnum_casing_panels::Int
: the number of panels between the leading and trailing edge of the duct on the duct inner side (casing)
Returns
vtan_tuple::NamedTuple
: a named tuple containing the body tangential surface velocities and various useful breakdowns thereof.
DuctAPE.get_body_tangential_velocities!
— Functionfunction getbodytangentialvelocities!( vtantuple, gamb, gamw, sigr, ivb, Vinf, totnode, totpanel, nnode, npanel, tangent, controlpoints, endpanelidxs, wakepanelidsalongcenterbodywakeinterface, wakepanelidsalongcasingwakeinterface, centerbodypanelidsalongcenterbodywakeinterface, ductpanelidsalongcasingwakeinterface, zpts, )
In-place version of get_body_tangential_velocities
.
Additional Arguments
zpts::NamedTuple
: a named tuple containing the z-coordinates of the control points of the duct casing, duct nacelle, and centerbody.
DuctAPE.calculate_vtheta
— Functioncalculate_vtheta(Gamma_tilde, r)
Calculate tangential velocity for a given net circulation and radial location
Arguments
Gamma_tilde::Matrix{Float}
: Sum of upstream circulation valuesr::Matrix{Float}
: blade element radial positions
DuctAPE.calculate_induced_velocities_on_bodywake
— Functioncalculate_induced_velocities_on_bodywake(
vz_w, vr_w, gamw, vz_r, vr_r, sigr, vz_b, vr_b, gamb, Vinf
)
Calculate the induced velocities on one of the body wakes (unit velocity inputs determine which one)
Arguments
vz_w::Matrix{Float}
: unit axial induced velocity of the wake onto the body wakevr_w::Matrix{Float}
: unit radial induced velocity of the wake onto the body wakegamw::Vector{Float}
: wake panel strengthsvz_r::Matrix{Float}
: unit axial induced velocity of the rotor onto the body wakevr_r::Matrix{Float}
: unit radial induced velocity of the rotor onto the body wakesigr::Vector{Float}
: rotor panel strengthsvz_b::Matrix{Float}
: unit axial induced velocity of the bodies onto the body wakevr_b::Matrix{Float}
: unit radial induced velocity of the bodies onto the body wakegamb::Vector{Float}
: body panel strengthsVinf::Vector{Float}
: one element vector containing the velocity magnitude
Pressures
DuctAPE.steady_cp
— Functionsteady_cp(Vs, Vinf, Vref)
Calculate steady pressure coefficients for a given surface velocity.
Arguments
Vs::Vector{Float}
: the surface velocitiesVinf::Vector{Float}
: one element vector with freestream mangnitudeVref::Vector{Float}
: one element vector with reference velocity used for non-dimensionalization
Returns
cp::Vector{Float}
: the steady pressure coefficients
DuctAPE.steady_cp!
— Functionsteady_cp!(cp, Vs, Vinf, Vref)
In-place verison of steady_cp
.
DuctAPE.calculate_entropy_jumps
— Functioncalculate_entropy_jumps(sigr, Cz_rotor)
Calculate jumps in entropy across the disks.
Arguments
sigr::Matrix{Float}
: rotor source panel strengthsCz_rotor::Vector{Float}
: absolute axial velocity on rotor blade elements
Returns
deltaS::Vector{Float}
: entropy jump across rotor disks
DuctAPE.calculate_rotor_jumps
— Functioncalculate_rotor_jumps(Gamr, Omega, B, sigr, Cz_rotor)
Calculate net circulation and enthalpy and entropy disk jumps
Arguments
Gamr::Matrix{Float}
: Blade element circulation strengthsOmega::Vector{Float}
: rotor rotation ratesB::Vector{Float}
: blade count for each rotor (usually integers but could be a float)sigr::Matrix{Float}
: rotor source panel strengthsCz_rotor::Vector{Float}
: absolute axial velocity on rotor blade elements
Returns
Gamma_tilde::Matrix{Float}
: net upstream circulationHtilde::Matrix{Float}
: net upstream enthalpy jumpsStilde::Matrix{Float}
: net upstream entropy jumps
DuctAPE.delta_cp
— Functiondelta_cp(deltaH, deltaS, Ctheta, Vref)
Calculate change in pressure coefficient aft of rotor, due to rotor
Arguments
deltaH::Vector{Float}
: Enthalpy jumps across disksdeltaS::Vector{Float}
: Entropy jumps across disks`Ctheta::Vector{Float}
: tangenetial velocityVref::Vector{Float}
: reference velocity for non-dimensionalization
Returns
delta_cp::Vector{Float}
: pressure rises due to rotor disks
DuctAPE.calculate_body_delta_cp!
— Functioncalculate_body_delta_cp!(cp, Gamr, sigr, Cz_rotor, Vref, Omega, B, cpr, casing_panel_ids_aft_of_rotors, centerbody_panel_ids_aft_of_rotors)
Augment surface pressure by change in pressure coefficient due to rotors specifically on the body panels aft of the rotors.
Arguments
cp::Vector{Float}
: steady pressure coeffients, modified in-place to include rotor effects.Gamr::Matrix{Float}
: Blade element circulation strengthssigr::Matrix{Float}
: rotor source panel strengthsCz_rotor::Vector{Float}
: absolute axial velocity on rotor blade elementsVref::Vector{Float}
: one element vector with reference velocity used for non-dimensionalizationOmega::Vector{Float}
: rotor rotation ratesB::Vector{Float}
: blade count for each rotor (usually integers but could be a float)cpr::Vector{Float}
: control point radial positions of body panelscasing_panel_ids_aft_of_rotors::Vector{Int}
: duct indices of control point radial positions aft of rotorscenterbody_panel_ids_aft_of_rotors::Vector{Int}
: centerbody indices of control point radial positions aft of rotors
DuctAPE.calculate_bodywake_delta_cp
— Functioncalculate_bodywake_delta_cp(Gamr, sigr, Cz_rotor, Vref, Omega, B, cpr; body="duct")
Calculate change in pressure coefficient due to rotors specifically on the body wakes
Arguments
Gamr::Matrix{Float}
: Blade element circulation strengthssigr::Matrix{Float}
: rotor source panel strengthsCz_rotor::Vector{Float}
: absolute axial velocity on rotor blade elementsVref::Vector{Float}
: one element vector with reference velocity used for non-dimensionalizationOmega::Vector{Float}
: rotor rotation ratesB::Vector{Float}
: blade count for each rotor (usually integers but could be a float)cpr::Vector{Float}
: control point radial positions of body wake "panels"
Keyword Arguments
body::String="duct"
: flag as to whether the body in question is a duct or centerbody.
DuctAPE.get_body_cps
— Functiongetbodycps( Vtanin, Vtanout, Gamr, sigr, Czrotor, Vinf, Vref, B, Omega, casingpanelidsaftofrotors, centerbodypanelidsaftof_rotors, controlpoints, endpanelidxs, zpts, )
Description
Arguments
Vtan_in::Vector{Float}
: Tangential velocity on the inside of the body panelsVtan_out::Vector{Float}
: Tangential velocity on the outside of the body panelsGamr::Matrix{Float}
: Blade element circulation strengthssigr::Matrix{Float}
: rotor source panel strengthsCz_rotor::Vector{Float}
: absolute axial velocity on rotor blade elementsVinf::Vector{Float}
: one element vector with freestream mangnitudeVref::Vector{Float}
: one element vector with reference velocity used for non-dimensionalizationB::Vector{Float}
: blade count for each rotor (usually integers but could be a float)Omega::Vector{Float}
: rotor rotation ratescasing_panel_ids_aft_of_rotors::Vector{Int}
: duct indices of control point radial positions aft of rotorscenterbody_panel_ids_aft_of_rotors::Vector{Int}
: centerbody indices of control point radial positions aft of rotorscontrolpoints::Matrix{Float}
: control point locations for each panelendpanelidxs::Matrix{Int}
: the indices of the first and last panels for each bodyzpts::NamedTuple
: a named tuple containing the z-coordinates of the control points of the duct casing, duct nacelle, and centerbody.
Returns
cp_tuple::NamedTuple
: body surface velocities and various useful breakdowns thereof.
DuctAPE.get_body_cps!
— Functionget_body_cps!(
cp_tuple,
Vtan_in,
Vtan_out,
Gamr,
sigr,
Cz_rotor,
Vinf,
Vref,
B,
Omega,
duct_panel_ids_aft_of_rotors,
centerbody_panel_ids_aft_of_rotors,
controlpoints,
endpanelidxs,
zpts,
)
In-place version of get_body_cps
.
DuctAPE.get_bodywake_cps
— Functionget_bodywake_cps(
Gamr,
vz_w,
vr_w,
gamw,
vz_r,
vr_r,
sigr,
vz_b,
vr_b,
gamb,
panels,
Cz_rotor,
Omega,
B,
Vinf,
Vref;
body="duct",
)
Calculate the pressure coefficient distributions on one of the body wakes
Arguments
Gamr::Matrix{Float}
: Blade element circulation strengthsvz_w::Matrix{Float}
: unit axial induced velocity of the wake onto the body wakevr_w::Matrix{Float}
: unit radial induced velocity of the wake onto the body wakegamw::Vector{Float}
: wake panel strengthsvz_r::Matrix{Float}
: unit axial induced velocity of the rotor onto the body wakevr_r::Matrix{Float}
: unit radial induced velocity of the rotor onto the body wakesigr::Vector{Float}
: rotor panel strengthsvz_b::Matrix{Float}
: unit axial induced velocity of the bodies onto the body wakevr_b::Matrix{Float}
: unit radial induced velocity of the bodies onto the body wakegamb::Vector{Float}
: body panel strengthspanels::NamedTuple
: A named tuple containing bodywake "panel" geometriesCz_rotor::Vector{Float}
: absolute axial velocity on rotor blade elementsOmega::Vector{Float}
: rotor rotation ratesB::Vector{Float}
: blade count for each rotor (usually integers but could be a float)Vinf::Vector{Float}
: one element vector containing the velocity magnitudeVref::Vector{Float}
: one element vector with reference velocity used for non-dimensionalization
Keyword Arguments
body::String="duct"
: flag as to whether the body in question is a duct or centerbody.
DuctAPE.forces_from_pressure
— Functionforces_from_pressure(cp_in, cp_out, panels; rhoinf=1.225, Vref=1.0)
Calculate dimensional and non-dimensional axial force on a single body
Arguments
cp_in::Vector{Float}
: pressure coefficient on inside of body surfacescp_out::Vector{Float}
: pressure coefficients on outside of body surfacespanels::NamedTuple
: A named tuple containing panel geometry information
Keyword Arguments
rhoinf::Float=1.225
: reference density for non-dimensionalizationVref::Float=1.0
: reference velocity for non-dimensionalization
Returns
thrust::Vector{Float}
: dimensional axial forceforce_coeff::Vector{Float}
: non-dimensional axial force
DuctAPE.forces_from_pressure!
— Functionforces_from_pressure!(CFx, cfx, cp_in, cp_out, panels; rhoinf=1.225, Vref=1.0)
In-place version of forces_from_pressure
.
DuctAPE.forces_from_TEpanels!
— Functionforces_from_TEpanels!(
thrust, force_coeff, cp_in, cp_out, panels; rhoinf=1.225, Vref=1.0
)
Add force induced by trailing edge gap panels to total forces.
Arguments
thrust::Vector{Float}
: dimensional axial forceforce_coeff::Vector{Float}
: non-dimensional axial forcecp_in::Vector{Float}
: pressure coefficient on inside of body surfacescp_out::Vector{Float}
: pressure coefficients on outside of body surfacespanels::NamedTuple
: A named tuple containing panel geometry information
Keyword Arguments
rhoinf::Float=1.225
: reference density for non-dimensionalizationVref::Float=1.0
: reference velocity for non-dimensionalization
Rotor Performance
DuctAPE.inviscid_rotor_thrust
— Functioninviscid_rotor_thrust(Ctheta_rotor, Gamma_tilde, rotor_panel_length, rhoinf)
Calculate inviscid rotor thrust.
Arguments
Ctheta_rotor::Vector{Float}
: Absolute tangential velocity on rotor blade elementsGamma_tilde::Matrix{Float}
: net upstream rotor circulationrotor_panel_length::Vector{Float}
: dimensional lengths on which blade element values applyrhoinf::Float
: freestream density
Returns
Tinv::Vector{Float}
: inviscid dimensional thrustdTi::Vector{Float}
: inviscid dimensional thrust distribution
DuctAPE.inviscid_rotor_thrust!
— Functioninviscid_rotor_thrust!(
Tinv, dTi, Ctheta_rotor, Gamma_tilde, rotor_panel_length, rhoinf
)
In-place version of inviscid_rotor_thrust
.
DuctAPE.viscous_rotor_thrust
— Functionviscous_rotor_thrust(
Cz_rotor, Cmag_rotor, B, chord, rotor_panel_length, cd, rhoinf
)
Calculate visous rotor "thrust."
Arguments
Cz_rotor::Vector{Float}
: Absolute axial velocity on rotor blade elementsCmag_rotor::Vector{Float}
: Absolute inflow velocity magnitude on rotor blade elementsB::Vector{Float}
: blade count for each rotor (usually integers but could be a float)chord::Vector{Float}
: blade element chord lengthsrotor_panel_length::Vector{Float}
: dimensional lengths on which blade element values applycd::Vector{Float}
: drag coefficient for each blade elementrhoinf::Float
: freestream density
Returns
Tvisc::Vector{Float}
: viscous dimensional thrustdTv::Vector{Float}
: viscous dimensional thrust distribution
DuctAPE.viscous_rotor_thrust!
— Functionviscous_rotor_thrust!(
Tvisc, dTv, Cz_rotor, Cmag_rotor, B, chord, rotor_panel_length, cd, rhoinf
)
In-place version of viscous_rotor_thrust
.
DuctAPE.inviscid_rotor_torque
— Functioninviscid_rotor_torque(
Cz_rotor, rotor_panel_center, rotor_panel_length, Gamma_tilde, rhoinf
)
Calculate inviscid rotor torque.
Arguments
Cz_rotor::Vector{Float}
: Absolute axial velocity on rotor blade elementsrotor_panel_center::Vector{Float}
: radial location of rotor blade elementsrotor_panel_length::Vector{Float}
: dimensional lengths on which blade element values applyGamma_tilde::Matrix{Float}
: net upstream rotor circulationrhoinf::Float
: freestream density
Returns
Qinv::Vector{Float}
: inviscid dimensional thrustdQi::Vector{Float}
: inviscid dimensional thrust distribution
DuctAPE.inviscid_rotor_torque!
— Functioninviscid_rotor_torque!(
Qinv, dQi, Cz_rotor, rotor_panel_center, rotor_panel_length, Gamma_tilde, rhoinf
)
In-place version of inviscid_rotor_torque
.
DuctAPE.viscous_rotor_torque
— Functionviscous_rotor_torque(
Ctheta_rotor, Cmag_rotor, B, chord, rotor_panel_center, rotor_panel_length, cd, rhoinf
)
Calculate viscous rotor torque.
Arguments
Ctheta_rotor::Vector{Float}
: Absolute tangential velocity on rotor blade elementsCmag_rotor::Vector{Float}
: Absolute inflow velocity magnitude on rotor blade elementsB::Vector{Float}
: blade count for each rotor (usually integers but could be a float)chord::Vector{Float}
: blade element chord lengthsrotor_panel_center::Vector{Float}
: radial location of rotor blade elementsrotor_panel_length::Vector{Float}
: dimensional lengths on which blade element values applycd::Vector{Float}
: drag coefficient for each blade elementrhoinf::Float
: freestream density
Returns
Qvisc::Vector{Float}
: viscous dimensional thrustdQv::Vector{Float}
: viscous dimensional thrust distribution
DuctAPE.viscous_rotor_torque!
— Functionviscous_rotor_torque!(
Qvisc,
dQv,
Ctheta_rotor,
Cmag_rotor,
B,
chord,
rotor_panel_center,
rotor_panel_length,
cd,
rhoinf
)
In-place version of viscous_rotor_torque
.
DuctAPE.rotor_power
— Functionrotor_power(Q, dQ, Omega)
Calculate power from torque and rotation rate.
Arguments
Q::Vector{Float}
: dimensional thrustdQ::Vector{Float}
: dimensional thrust distributionOmega::Vector{Float}
: rotor rotation rates
Returns
P::Vector{Float}
: dimensional powerdP::Vector{Float}
: dimensional thrust distribution
DuctAPE.rotor_power!
— Functionrotor_power!(P, dP, Q, dQ, Omega)
In-place version of rotor_power
.
DuctAPE.get_total_efficiency
— Functionget_total_efficiency(total_thrust, total_power, Vinf)
Get total efficiency.
Arguments
total_thrust::Vector{Float}
: total thrusttotal_power::Vector{Float}
: total powerVinf::Vector{Float}
: one element vector freestream velocity magnitude
Returns
- `total_efficiency::Vector{Float} : total efficiency
DuctAPE.get_total_efficiency!
— Functionget_total_efficiency!(eta, total_thrust, total_power, Vinf)
In-place version of get_total_efficiency
.
DuctAPE.get_induced_efficiency
— Functionget_induced_efficiency(Tinv, Tduct, Pinv, Vinf)
Get rotor efficiency induced by presence of the duct.
Arguments
Tinv::Vector{Float}
: inviscid dimensional thrustTduct::Vector{Float}
: duct thrustPinv::Vector{Float}
: inviscid dimensional powerVinf::Vector{Float}
: one element vector freestream velocity magnitude
Returns
induced_efficiency::Vector{Float}
: rotor efficiency induced by duct
DuctAPE.get_induced_efficiency!
— Functionget_induced_efficiency!(eta_inv, Tinv, Tduct, Pinv, Vinf)
In-place version of get_induced_efficiency
.
DuctAPE.get_ideal_efficiency
— Functionget_ideal_efficiency(total_thrust, rhoinf, Vinf, Rref)
Compute ducted fan ideal efficiency
Arguments
total_thrust::Vector{Float}
: total thrust from rotors and ductrhoinf::Float
: freestream densityVinf::Vector{Float}
: one element vector freestream velocity magnitudeRref::Vector{Float}
: one element vector reference rotor tip radius
Returns
ideal_efficiency::Vector{Float}
: ideal ducted fan efficiency
DuctAPE.tqpcoeff
— Functiontqpcoeff(thrust, torque, power, rhoinf, Omega, Rref)
Calculate non-dimensional thrust, torque, and power coefficients
Arguments
thrust::Vector{Float}
: dimensional thrusttorque::Vector{Float}
: dimensional torquepower::Vector{Float}
: dimensional powerrhoinf::Float
: freestream densityOmega::Vector{Float}
: rotor rotation ratesRref::Vector{Float}
: one element vector reference rotor tip radius
Returns
CT::Vector{Float}
: thrust coefficientCQ::Vector{Float}
: torque coefficientCP::Vector{Float}
: power coefficient
DuctAPE.tqpcoeff!
— Functiontqpcoeff!(CT, CQ, CP, thrust, torque, power, rhoinf, Omega, Rref)
In-place version of tqpcoeff
.
DuctAPE.get_blade_loads
— Functionget_blade_loads(Cmag_rotor, beta1, cl, cd, chords, rhoinf)
Get loading along blades.
Arguments
Cmag_rotor::Vector{Float}
: blade element inflow magnitudesbeta1::Vector{Float}
: blade element inflow anglescl::Vector{Float}
: blade element lift coefficientscd::Vector{Float}
: blade element drag coefficientschords::Vector{Float}
: blade element chord lengthsrhoinf::Vector{Float}
: one element freestream density
Returns
Np::Vector{Float}
: blade loading per unit length in the normal direction: N'Tp::Vector{Float}
: blade loading per unit length in the tangential direction: T'
DuctAPE.get_blade_loads!
— Functionget_blade_loads!(Np, Tp, Cmag_rotor, beta1, cl, cd, chords, rhoinf, cache)
In-place version of get_blade_loads
.
Boundary Layer
Thermodynamics
DuctAPE.sa1
— Functionsa1(altitude; hardness=50)
Standard atmosphere temperature and pressure in SI units blended between the first linear portion and the constant portion.
DuctAPE.sa2
— Functionsa2(altitude; hardness=50)
Standard atmosphere temperature and pressure in SI units blended between the the constant portion and the second linear portion.
DuctAPE.standard_atmosphere
— Functionstandard_atmosphere(altitude; hardness=25)
standard_atmosphere(::Imperial, altitude; hardness=25)
Smoothed fits to the Standard Atmosphere model.
Assumes calorically imperfect gas.
Arguments
altitude::Float
: Altitude in meters for SI units, or feet for Imperial units
Keyword Arguments:
hardness::float
: hardness factor for sigmoid blend
Returns
static_temperature::Float
: Static temperaturestatic_pressure::Float
: Static pressurestatic_density::Float
: Static densitystatic_dynamic_viscosity::Float
: Static dynamic Viscosity
DuctAPE.ideal_gas_rho
— Functionideal_gas_rho(P, T)
Ideal gas law for calculating density in SI units.
DuctAPE.sutherlands_law
— Functionsutherlands_law(
static_temperure, mu_sea_level=1.789e-5, T_sea_level=288.15, S=110.4
)
Sutherland's law in SI units for calculating air viscosity relative to sea level.
DuctAPE.speed_of_sound
— Functionspeed_of_sound(static_pressure, static_density; gamma=1.4)
Speed of sound from isentropic relations
DuctAPE.calculate_mach
— Functioncalculate_mach(edge_velocity, speed_of_sound)
Mach number from velocity and speed of sound
DuctAPE.total_temperature
— Functiontotal_temperature(static_temperature, M; gamma=1.4)
Total temperature from isentropic relations
DuctAPE.total_pressure
— Functiontotal_pressure(static_pressure, M; gamma=1.4)
Total pressure from isentropic relations
DuctAPE.static_temperature
— Functionstatic_temperature(total_temperature, M; gamma=1.4)
Static temperature from isentropic relations
DuctAPE.static_pressure
— Functionstatic_pressure(total_pressure, M; gamma=1.4)
Static pressure from isentropic relations
DuctAPE.static_density
— Functionstatic_density(static_pressure, speed_of_sound; gamma=1.4)
Static density from isentropic relations
DuctAPE.convert_temperature_to_kelvin
— Functionconvert_temperature_to_kelvin(::Units, T)
Convert from Fahrenheit to Kelvin or Return temperature if already in SI units
DuctAPE.convert_viscosity
— Functionconvert_viscosity(::SI, mu)
Convert viscosity from Imperial units to SI or Return input if already SI units.
General Boundary Layer Functions
DuctAPE.arc_lengths_from_panel_lengths
— Functionarc_lengths_from_panel_lengths(duct_panel_lengths, bl_ids)
Cumulative sum of panel lengths for the given section of surface associated with the upper or lower boundary layer.
Arguments:
duct_panel_lengths::Vector{Float}
: vector of panel lengths (called influencelength in bodyvortex_panels) associated with the duct (casing + nacelle).
Returns:
s::Vector{Float}
: cumulative sum of panel lengths between control points in the given index range, starting from zero.
DuctAPE.split_at_stagnation_point
— Functionsplitatstagnationpoint(ductpanellengths, npanels_casing)
Split the duct body surface at the leading edge of the duct.
Arguments:
duct_panel_lengths::Vector{Float}
: Vector of panel lengths for the duct from casing trailing edge clockwise to nacelle trailing edge.n_panels_casing::Int
: number of panels comprising the casing side of the duct
Returns:
s_upper::Vector{Float}
: cumulative sum of upper side (nacelle) panel lengthss_lower::Vector{Float}
: cumulative sum of lower side (casing) panel lengths
DuctAPE.bl_step_fun
— Functionbl_step_fun(n, m, p)
Function used in determining step sizes for boundary layer calculation. f(n) = m*n^p
Given a number of steps, n ∈ [1:N], provides the cumulative step lengths according to the power, p, and the multiplicative factor, m; where p determined from the set_boundary_layer_steps
functions.
DuctAPE.set_boundary_layer_steps
— Functionset_boundary_layer_steps(N::Int, first_step_size, total_length)
Sets boundary layer steps based on desired number of steps (must be an Integer), an initial step size, and the total cumulative length of the steps.
Arguments:
N::Int
: Number of steps to takefirst_step_size::Float
: size of first step (which ism
inbl_step_fun
)total_length::Float
: total surface length to divide up.
Returns:
steps::Vector{Float}
: steps along surface length satisfying the equation: f(n) = m*n^p with the condition thatm
is the first step size and f(N) =total_length
DuctAPE.RK2
— FunctionRK2(f, y, s, ds, parameters)
2nd Order Runge-Kutta integration scheme.
Arguments:
f::function_handle
: residual function for integrationy::Vector{Float}
: current statess::Float
: current positionds::Float
: step sizeparameters::NamedTuple
: BoundaryLayerOptions and other various parameters
DuctAPE.RK4
— FunctionRK4(f, y, s, ds, parameters)
4th Order Runge-Kutta integration scheme.
Arguments:
f::function_handle
: residual function for integrationy::Vector{Float}
: current statess::Float
: current positionds::Float
: step sizeparameters::NamedTuple
: BoundaryLayerOptions and other various parameters
Head's Method Specific Functions
DuctAPE.setup_boundary_layer_functions_head
— Functionsetup_boundary_layer_functions_head(
s,
vtan_duct,
duct_control_points,
operating_point,
boundary_layer_options;
verbose=false
)
Arguments:
s::Vector{Float}
: cumulative sum of panel lengths between control points in the given index range, starting from zero.vtan_duct::Vector{Float}
: tangential velocity magnitudes for the entire ductduct_control_points::Matrix{Float}
: Control point coordinates along the duct surfaceoperating_point::OperatingPoint
: OperatingPoint objectboundary_layer_options::BoundaryLayerOptions
: BoundaryLayerOptions object
Returns:
boundary_layer_parameters::NamedTuple
: Namped Tuple containing boundary layer solver parameters:edge_velocity::FLOWMath.Akima
: spline of edge velocities relative to surface lengthedge_acceleration::FLOWMath.Akima
: spline of edge acceleration (dUe/ds) relative to surface lengthedge_density::FLOWMath.Akima
: spline of edge density relative to surface lengthedge_viscosity::FLOWMath.Akima
: spline of edge viscosity relative to surface length
DuctAPE.calculate_H
— Functioncalculate_H(H1)
Calculate the value of the shape factor used in Head's method.
DuctAPE.calculate_cf
— Functioncalculate_cf(H, Red2)
Calculate the skin friction coefficient used in Head's method
DuctAPE.boundary_layer_residual_head
— Functionboundary_layer_residual_head(y, parameters, s)
Out of place residual function for Head's method.
DuctAPE.boundary_layer_residual_head!
— Functionboundary_layer_residual_head!(dy, y, parameters, s)
Calculate dy give the current states, y, the input position, s, and various parameters.
DuctAPE.solve_head_boundary_layer!
— Functionsolve_head_boundary_layer!(f, ode, initial_states, steps, parameters; verbose=false)
Integrate the turbulent boundary layer using a Runge-Kutta method.
Arguments:
f::function_handle
: Governing residual equations to integrateode::function_handle
: ODE method to use (RK2 or RK4)initial_states::Float
: initial statessteps::Vector{Float}
: steps for integrationparameters::NamedTuple
: boundary layer solve options and other parameters
solve_head_boundary_layer!(::DiffEq, ode, initial_states, steps, parameters; verbose=false)
Integrate the turbulent boundary layer using a Runge-Kutta method.
Arguments:
f::function_handle
: Governing residual equations to integrateode::function_handle
: ODE method to use (one of the DifferentialEquations.jl options)initial_states::Float
: initial statessteps::Vector{Float}
: steps for integrationparameters::NamedTuple
: boundary layer solve options and other parameters
Viscous Drag
DuctAPE.squire_young
— Functionsquire_young(d2, Ue, Uinf, H12)
Squire-Young formula for the viscous drag coeffiecient of one side of a 2D body.
Arguments:
d2::Float
: Momentum thickness at separation extrapolated back to the trailing edge (d2 = d2sep+rsep-rTE)Ue::Float
: Edge velocity at separation pointUinf::Float
: Freestream velocityH12::Float
: Boundary layer shape factor at separation point
Note: if no separation occurs, the inputs are simply the final values for the boundary layer.
Returns:
cdc::Float
: viscous drag coefficient times reference chord
DuctAPE.total_viscous_drag_duct
— Functiontotal_viscous_drag_duct(cd_upper, cd_lower, rotor_tip_radius, Vref, rhoinf)
Calculate the total viscous drag of the duct from Squire-Young drag coefficients, integrated about exit circumference.
Arguments:
cdc_upper::Float
: upper side drag coefficient times refernce chordcdc_lower::Float
: lower side drag coefficient times refernce chordrotor_tip_radius::Float
: radius used for integrating circumferentiallyVref::Float
: reference velocity (Vinf)rhoinf::Float
: freestream density
Returns:
viscous_drag::Float
: viscous drag on duct
DuctAPE.compute_viscous_drag_duct
— Functioncompute_viscous_drag_duct(
boundary_layer_options::BoundaryLayerOptions,
Vtan_duct,
cp_duct,
duct_panel_lengths,
rotor_tip_radius,
operating_point,
reference_parameters,
)
Determine total, dimensional viscous drag on the duct.
Arguments:
boundary_layer_options::BoundaryLayerOptions
: BoundaryLayerOptions object, used for dispatch as wellVtan_duct::Vector{Float}
: tangential velocity magnitudes for the entire ductduct_panel_lengths::Vector{Float}
: panel lengths for the entire ductrotor_tip_radius::Float
: radius at duct trailing edge (casing side)operating_point::OperatingPoint
: OperatingPoint objectreference_parameters::ReferenceParameters
: ReferenceParameters object
Returns:
duct_viscous_drag::Float
: total viscous drag of ductboundary_layer_outputs::NamedTuple
: named tuple of various boundary layer related outputs:stagnation_indices
upper_solved_states
upper_solved_steps
lower_solved_states
lower_solved_steps
surface_length_upper
surface_length_lower
split_ratio
separation_point_ratio_upper
separation_point_ratio_lower
cdc_upper
cdc_lower
DuctAPE.compute_single_side_drag_coefficient
— Functioncompute_single_side_drag_coefficient(
steps,
rotor_tip_radius,
operating_point,
boundary_layer_options;
verbose=false,
)
Solve integral boundary layer and obtain viscous drag coefficient from Squire-Young formula for one side of the duct (side being defined as portion of surface on once side of the stagnation point)
Arguments:
steps::Vector{Float}
: positions along surface for integrationrotor_tip_radius::Float
: radius at duct trailing edge (casing side)Vref::Float
: reference velocitysetup_boundary_layer_functions::NamedTuple
: Various Akima splines and other functions for boundary layer valuesboundary_layer_options::BoundaryLayerOptions
: BoundaryLayerOptions object
Returns:
cd::Float
: viscous drag coefficient
DuctAPE.compute_viscous_drag_duct_schlichting
— Functioncompute_viscous_drag_duct_schlichting(
vtan_duct_TE, duct_chord, TE_radius, operating_point
)
Computes Schlichting approximation of skin friction drag dimensionalized to drag per unit length using duct chord, and using the trailing edge circuference as the total length.
Arguments:
vtan_duct_TE::Float
: tangential velocity at the duct trailing edgeduct_chord::Float
: length of ductTE_radius::Float
: radius of the trailing edge pointoperating_point::OperatingPoint
: OperatingPoint object