DuctAPE.post_processFunction
post_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 variables
  • prepost_containers::NamedTuple : the named tuple containing pre-allocated containers for the pre- and post-processing intermediate calculations
  • solve_container_cache::NamedTuple : the cache and dimensions for intermediate values in the residual calculation
  • solve_parameter_cache_vector::Vector{Float} : the applicably typed cache vector for the solve parameters
  • solve_parameter_cache_dims::NamedTuple : the dimensions of the solver parameters
  • operating_point::OperatingPoint : the operating point being analyzed
  • reference_parameters::ReferenceParameters : a ReferenceParameters object
  • BoundaryLayerOptions::BoundaryLayerOptions : a BoundaryLayerOptions object
  • A_bb_LU::LinearAlgebra.LU : LinearAlgebra LU factorization of the LHS matrix
  • airfoils::Vector{AFType} : A matrix of airfoil types associated with each of the blade elements
  • idmaps::NamedTuple : A named tuple containing index mapping used in bookkeeping throughout solve and post-process
  • problem_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 written
  • checkoutfileexists=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_uppercdc_upper`
      • cdc_lowercdc_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
source

Velocities

DuctAPE.get_body_tangential_velocitiesFunction
get_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 strengths
  • gamw::Vector{Float} : the wake panel strengths
  • sigr::Vector{Float} : the rotor panel strengths
  • ivb::NamedTuple : the unit induced velocities on the bodies
  • Vinf::Vector{Float} : one element vector containing the freestream magnitude
  • totnode::Int : total number of nodes between all bodies
  • totpanel::Int : total number of panels between all bodies
  • nnode::Vector{Int} : number of nodes in each body
  • npanel::Vector{Int} : number of panels in each body.
  • tangent::Matrix{Float} : unit tangent vectors for each panel
  • controlpoints::Matrix{Float} : control point locations for each panel
  • endpanelidxs::Matrix{Int} : the indices of the first and last panels for each body
  • wake_panel_ids_along_centerbody_wake_interface::Vector{Int} : the indices of the wake panels coincident with the centerbody panels
  • wake_panel_ids_along_casing_wake_interface::Vector{Int} : the indices of the wake panels coincident with the duct casing (inner surface) panels
  • centerbody_panel_ids_along_centerbody_wake_interface::Vector{Int} : the indices of the centerbody panels coincident with the wake panels
  • duct_panel_ids_along_casing_wake_interface::Vector{Int} : the indices of the duct panels coincident with the wake panels
  • num_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.
source
DuctAPE.get_body_tangential_velocities!Function

function 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.
source
DuctAPE.calculate_vthetaFunction
calculate_vtheta(Gamma_tilde, r)

Calculate tangential velocity for a given net circulation and radial location

Arguments

  • Gamma_tilde::Matrix{Float} : Sum of upstream circulation values
  • r::Matrix{Float} : blade element radial positions
source
DuctAPE.calculate_induced_velocities_on_bodywakeFunction
calculate_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 wake
  • vr_w::Matrix{Float} : unit radial induced velocity of the wake onto the body wake
  • gamw::Vector{Float} : wake panel strengths
  • vz_r::Matrix{Float} : unit axial induced velocity of the rotor onto the body wake
  • vr_r::Matrix{Float} : unit radial induced velocity of the rotor onto the body wake
  • sigr::Vector{Float} : rotor panel strengths
  • vz_b::Matrix{Float} : unit axial induced velocity of the bodies onto the body wake
  • vr_b::Matrix{Float} : unit radial induced velocity of the bodies onto the body wake
  • gamb::Vector{Float} : body panel strengths
  • Vinf::Vector{Float} : one element vector containing the velocity magnitude
source

Pressures

DuctAPE.steady_cpFunction
steady_cp(Vs, Vinf, Vref)

Calculate steady pressure coefficients for a given surface velocity.

Arguments

  • Vs::Vector{Float} : the surface velocities
  • Vinf::Vector{Float} : one element vector with freestream mangnitude
  • Vref::Vector{Float} : one element vector with reference velocity used for non-dimensionalization

Returns

  • cp::Vector{Float} : the steady pressure coefficients
source
DuctAPE.calculate_entropy_jumpsFunction
calculate_entropy_jumps(sigr, Cz_rotor)

Calculate jumps in entropy across the disks.

Arguments

  • sigr::Matrix{Float} : rotor source panel strengths
  • Cz_rotor::Vector{Float} : absolute axial velocity on rotor blade elements

Returns

  • deltaS::Vector{Float} : entropy jump across rotor disks
source
DuctAPE.calculate_rotor_jumpsFunction
calculate_rotor_jumps(Gamr, Omega, B, sigr, Cz_rotor)

Calculate net circulation and enthalpy and entropy disk jumps

Arguments

  • Gamr::Matrix{Float} : Blade element circulation strengths
  • Omega::Vector{Float} : rotor rotation rates
  • B::Vector{Float} : blade count for each rotor (usually integers but could be a float)
  • sigr::Matrix{Float} : rotor source panel strengths
  • Cz_rotor::Vector{Float} : absolute axial velocity on rotor blade elements

Returns

  • Gamma_tilde::Matrix{Float} : net upstream circulation
  • Htilde::Matrix{Float} : net upstream enthalpy jumps
  • Stilde::Matrix{Float} : net upstream entropy jumps
source
DuctAPE.delta_cpFunction
delta_cp(deltaH, deltaS, Ctheta, Vref)

Calculate change in pressure coefficient aft of rotor, due to rotor

Arguments

  • deltaH::Vector{Float} : Enthalpy jumps across disks
  • deltaS::Vector{Float} : Entropy jumps across disks`
  • Ctheta::Vector{Float} : tangenetial velocity
  • Vref::Vector{Float} : reference velocity for non-dimensionalization

Returns

  • delta_cp::Vector{Float} : pressure rises due to rotor disks
source
DuctAPE.calculate_body_delta_cp!Function
calculate_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 strengths
  • sigr::Matrix{Float} : rotor source panel strengths
  • Cz_rotor::Vector{Float} : absolute axial velocity on rotor blade elements
  • Vref::Vector{Float} : one element vector with reference velocity used for non-dimensionalization
  • Omega::Vector{Float} : rotor rotation rates
  • B::Vector{Float} : blade count for each rotor (usually integers but could be a float)
  • cpr::Vector{Float} : control point radial positions of body panels
  • casing_panel_ids_aft_of_rotors::Vector{Int} : duct indices of control point radial positions aft of rotors
  • centerbody_panel_ids_aft_of_rotors::Vector{Int} : centerbody indices of control point radial positions aft of rotors
source
DuctAPE.calculate_bodywake_delta_cpFunction
calculate_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 strengths
  • sigr::Matrix{Float} : rotor source panel strengths
  • Cz_rotor::Vector{Float} : absolute axial velocity on rotor blade elements
  • Vref::Vector{Float} : one element vector with reference velocity used for non-dimensionalization
  • Omega::Vector{Float} : rotor rotation rates
  • B::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.
source
DuctAPE.get_body_cpsFunction

getbodycps( 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 panels
  • Vtan_out::Vector{Float} : Tangential velocity on the outside of the body panels
  • Gamr::Matrix{Float} : Blade element circulation strengths
  • sigr::Matrix{Float} : rotor source panel strengths
  • Cz_rotor::Vector{Float} : absolute axial velocity on rotor blade elements
  • Vinf::Vector{Float} : one element vector with freestream mangnitude
  • Vref::Vector{Float} : one element vector with reference velocity used for non-dimensionalization
  • B::Vector{Float} : blade count for each rotor (usually integers but could be a float)
  • Omega::Vector{Float} : rotor rotation rates
  • casing_panel_ids_aft_of_rotors::Vector{Int} : duct indices of control point radial positions aft of rotors
  • centerbody_panel_ids_aft_of_rotors::Vector{Int} : centerbody indices of control point radial positions aft of rotors
  • controlpoints::Matrix{Float} : control point locations for each panel
  • endpanelidxs::Matrix{Int} : the indices of the first and last panels for each body
  • zpts::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.
source
DuctAPE.get_body_cps!Function
get_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.

source
DuctAPE.get_bodywake_cpsFunction
get_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 strengths
  • vz_w::Matrix{Float} : unit axial induced velocity of the wake onto the body wake
  • vr_w::Matrix{Float} : unit radial induced velocity of the wake onto the body wake
  • gamw::Vector{Float} : wake panel strengths
  • vz_r::Matrix{Float} : unit axial induced velocity of the rotor onto the body wake
  • vr_r::Matrix{Float} : unit radial induced velocity of the rotor onto the body wake
  • sigr::Vector{Float} : rotor panel strengths
  • vz_b::Matrix{Float} : unit axial induced velocity of the bodies onto the body wake
  • vr_b::Matrix{Float} : unit radial induced velocity of the bodies onto the body wake
  • gamb::Vector{Float} : body panel strengths
  • panels::NamedTuple : A named tuple containing bodywake "panel" geometries
  • Cz_rotor::Vector{Float} : absolute axial velocity on rotor blade elements
  • Omega::Vector{Float} : rotor rotation rates
  • B::Vector{Float} : blade count for each rotor (usually integers but could be a float)
  • Vinf::Vector{Float} : one element vector containing the velocity magnitude
  • Vref::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.
source
DuctAPE.forces_from_pressureFunction
forces_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 surfaces
  • cp_out::Vector{Float} : pressure coefficients on outside of body surfaces
  • panels::NamedTuple : A named tuple containing panel geometry information

Keyword Arguments

  • rhoinf::Float=1.225 : reference density for non-dimensionalization
  • Vref::Float=1.0 : reference velocity for non-dimensionalization

Returns

  • thrust::Vector{Float} : dimensional axial force
  • force_coeff::Vector{Float} : non-dimensional axial force
source
DuctAPE.forces_from_TEpanels!Function
forces_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 force
  • force_coeff::Vector{Float} : non-dimensional axial force
  • cp_in::Vector{Float} : pressure coefficient on inside of body surfaces
  • cp_out::Vector{Float} : pressure coefficients on outside of body surfaces
  • panels::NamedTuple : A named tuple containing panel geometry information

Keyword Arguments

  • rhoinf::Float=1.225 : reference density for non-dimensionalization
  • Vref::Float=1.0 : reference velocity for non-dimensionalization
source

Rotor Performance

DuctAPE.inviscid_rotor_thrustFunction
inviscid_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 elements
  • Gamma_tilde::Matrix{Float} : net upstream rotor circulation
  • rotor_panel_length::Vector{Float} : dimensional lengths on which blade element values apply
  • rhoinf::Float : freestream density

Returns

  • Tinv::Vector{Float} : inviscid dimensional thrust
  • dTi::Vector{Float} : inviscid dimensional thrust distribution
source
DuctAPE.inviscid_rotor_thrust!Function
inviscid_rotor_thrust!(
    Tinv, dTi, Ctheta_rotor, Gamma_tilde, rotor_panel_length, rhoinf
)

In-place version of inviscid_rotor_thrust.

source
DuctAPE.viscous_rotor_thrustFunction
viscous_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 elements
  • Cmag_rotor::Vector{Float} : Absolute inflow velocity magnitude on rotor blade elements
  • B::Vector{Float} : blade count for each rotor (usually integers but could be a float)
  • chord::Vector{Float} : blade element chord lengths
  • rotor_panel_length::Vector{Float} : dimensional lengths on which blade element values apply
  • cd::Vector{Float} : drag coefficient for each blade element
  • rhoinf::Float : freestream density

Returns

  • Tvisc::Vector{Float} : viscous dimensional thrust
  • dTv::Vector{Float} : viscous dimensional thrust distribution
source
DuctAPE.viscous_rotor_thrust!Function
viscous_rotor_thrust!(
    Tvisc, dTv, Cz_rotor, Cmag_rotor, B, chord, rotor_panel_length, cd, rhoinf
)

In-place version of viscous_rotor_thrust.

source
DuctAPE.inviscid_rotor_torqueFunction
inviscid_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 elements
  • rotor_panel_center::Vector{Float} : radial location of rotor blade elements
  • rotor_panel_length::Vector{Float} : dimensional lengths on which blade element values apply
  • Gamma_tilde::Matrix{Float} : net upstream rotor circulation
  • rhoinf::Float : freestream density

Returns

  • Qinv::Vector{Float} : inviscid dimensional thrust
  • dQi::Vector{Float} : inviscid dimensional thrust distribution
source
DuctAPE.inviscid_rotor_torque!Function
inviscid_rotor_torque!(
    Qinv, dQi, Cz_rotor, rotor_panel_center, rotor_panel_length, Gamma_tilde, rhoinf
)

In-place version of inviscid_rotor_torque.

source
DuctAPE.viscous_rotor_torqueFunction
viscous_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 elements
  • Cmag_rotor::Vector{Float} : Absolute inflow velocity magnitude on rotor blade elements
  • B::Vector{Float} : blade count for each rotor (usually integers but could be a float)
  • chord::Vector{Float} : blade element chord lengths
  • rotor_panel_center::Vector{Float} : radial location of rotor blade elements
  • rotor_panel_length::Vector{Float} : dimensional lengths on which blade element values apply
  • cd::Vector{Float} : drag coefficient for each blade element
  • rhoinf::Float : freestream density

Returns

  • Qvisc::Vector{Float} : viscous dimensional thrust
  • dQv::Vector{Float} : viscous dimensional thrust distribution
source
DuctAPE.viscous_rotor_torque!Function
viscous_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.

source
DuctAPE.rotor_powerFunction
rotor_power(Q, dQ, Omega)

Calculate power from torque and rotation rate.

Arguments

  • Q::Vector{Float} : dimensional thrust
  • dQ::Vector{Float} : dimensional thrust distribution
  • Omega::Vector{Float} : rotor rotation rates

Returns

  • P::Vector{Float} : dimensional power
  • dP::Vector{Float} : dimensional thrust distribution
source
DuctAPE.get_total_efficiencyFunction
get_total_efficiency(total_thrust, total_power, Vinf)

Get total efficiency.

Arguments

  • total_thrust::Vector{Float} : total thrust
  • total_power::Vector{Float} : total power
  • Vinf::Vector{Float} : one element vector freestream velocity magnitude

Returns

  • `total_efficiency::Vector{Float} : total efficiency
source
DuctAPE.get_induced_efficiencyFunction
get_induced_efficiency(Tinv, Tduct, Pinv, Vinf)

Get rotor efficiency induced by presence of the duct.

Arguments

  • Tinv::Vector{Float} : inviscid dimensional thrust
  • Tduct::Vector{Float} : duct thrust
  • Pinv::Vector{Float} : inviscid dimensional power
  • Vinf::Vector{Float} : one element vector freestream velocity magnitude

Returns

  • induced_efficiency::Vector{Float} : rotor efficiency induced by duct
source
DuctAPE.get_ideal_efficiencyFunction
get_ideal_efficiency(total_thrust, rhoinf, Vinf, Rref)

Compute ducted fan ideal efficiency

Arguments

  • total_thrust::Vector{Float} : total thrust from rotors and duct
  • rhoinf::Float : freestream density
  • Vinf::Vector{Float} : one element vector freestream velocity magnitude
  • Rref::Vector{Float} : one element vector reference rotor tip radius

Returns

  • ideal_efficiency::Vector{Float} : ideal ducted fan efficiency
source
DuctAPE.tqpcoeffFunction
tqpcoeff(thrust, torque, power, rhoinf, Omega, Rref)

Calculate non-dimensional thrust, torque, and power coefficients

Arguments

  • thrust::Vector{Float} : dimensional thrust
  • torque::Vector{Float} : dimensional torque
  • power::Vector{Float} : dimensional power
  • rhoinf::Float : freestream density
  • Omega::Vector{Float} : rotor rotation rates
  • Rref::Vector{Float} : one element vector reference rotor tip radius

Returns

  • CT::Vector{Float} : thrust coefficient
  • CQ::Vector{Float} : torque coefficient
  • CP::Vector{Float} : power coefficient
source
DuctAPE.tqpcoeff!Function
tqpcoeff!(CT, CQ, CP, thrust, torque, power, rhoinf, Omega, Rref)

In-place version of tqpcoeff.

source
DuctAPE.get_blade_loadsFunction
get_blade_loads(Cmag_rotor, beta1, cl, cd, chords, rhoinf)

Get loading along blades.

Arguments

  • Cmag_rotor::Vector{Float} : blade element inflow magnitudes
  • beta1::Vector{Float} : blade element inflow angles
  • cl::Vector{Float} : blade element lift coefficients
  • cd::Vector{Float} : blade element drag coefficients
  • chords::Vector{Float} : blade element chord lengths
  • rhoinf::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'
source
DuctAPE.get_blade_loads!Function
get_blade_loads!(Np, Tp, Cmag_rotor, beta1, cl, cd, chords, rhoinf, cache)

In-place version of get_blade_loads.

source

Boundary Layer

Thermodynamics

DuctAPE.sa1Function
sa1(altitude; hardness=50)

Standard atmosphere temperature and pressure in SI units blended between the first linear portion and the constant portion.

source
DuctAPE.sa2Function
sa2(altitude; hardness=50)

Standard atmosphere temperature and pressure in SI units blended between the the constant portion and the second linear portion.

source
DuctAPE.standard_atmosphereFunction
standard_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 temperature
  • static_pressure::Float : Static pressure
  • static_density::Float : Static density
  • static_dynamic_viscosity::Float : Static dynamic Viscosity
source
DuctAPE.sutherlands_lawFunction
sutherlands_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.

source

General Boundary Layer Functions

DuctAPE.arc_lengths_from_panel_lengthsFunction
arc_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.
source
DuctAPE.split_at_stagnation_pointFunction

splitatstagnationpoint(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 lengths
  • s_lower::Vector{Float} : cumulative sum of lower side (casing) panel lengths
source
DuctAPE.bl_step_funFunction
bl_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.

source
DuctAPE.set_boundary_layer_stepsFunction
set_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 take
  • first_step_size::Float : size of first step (which is m in bl_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 that m is the first step size and f(N) = total_length
source
DuctAPE.RK2Function
RK2(f, y, s, ds, parameters)

2nd Order Runge-Kutta integration scheme.

Arguments:

  • f::function_handle : residual function for integration
  • y::Vector{Float} : current states
  • s::Float : current position
  • ds::Float : step size
  • parameters::NamedTuple : BoundaryLayerOptions and other various parameters
source
DuctAPE.RK4Function
RK4(f, y, s, ds, parameters)

4th Order Runge-Kutta integration scheme.

Arguments:

  • f::function_handle : residual function for integration
  • y::Vector{Float} : current states
  • s::Float : current position
  • ds::Float : step size
  • parameters::NamedTuple : BoundaryLayerOptions and other various parameters
source

Head's Method Specific Functions

DuctAPE.setup_boundary_layer_functions_headFunction
setup_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 duct
  • duct_control_points::Matrix{Float} : Control point coordinates along the duct surface
  • operating_point::OperatingPoint : OperatingPoint object
  • boundary_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 length
    • edge_acceleration::FLOWMath.Akima : spline of edge acceleration (dUe/ds) relative to surface length
    • edge_density::FLOWMath.Akima : spline of edge density relative to surface length
    • edge_viscosity::FLOWMath.Akima : spline of edge viscosity relative to surface length
source
DuctAPE.solve_head_boundary_layer!Function
solve_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 integrate
  • ode::function_handle : ODE method to use (RK2 or RK4)
  • initial_states::Float : initial states
  • steps::Vector{Float} : steps for integration
  • parameters::NamedTuple : boundary layer solve options and other parameters
source
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 integrate
  • ode::function_handle : ODE method to use (one of the DifferentialEquations.jl options)
  • initial_states::Float : initial states
  • steps::Vector{Float} : steps for integration
  • parameters::NamedTuple : boundary layer solve options and other parameters
source

Viscous Drag

DuctAPE.squire_youngFunction
squire_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 point
  • Uinf::Float : Freestream velocity
  • H12::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
source
DuctAPE.total_viscous_drag_ductFunction
total_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 chord
  • cdc_lower::Float : lower side drag coefficient times refernce chord
  • rotor_tip_radius::Float : radius used for integrating circumferentially
  • Vref::Float : reference velocity (Vinf)
  • rhoinf::Float : freestream density

Returns:

  • viscous_drag::Float : viscous drag on duct
source
DuctAPE.compute_viscous_drag_ductFunction
compute_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 well
  • Vtan_duct::Vector{Float} : tangential velocity magnitudes for the entire duct
  • duct_panel_lengths::Vector{Float} : panel lengths for the entire duct
  • rotor_tip_radius::Float : radius at duct trailing edge (casing side)
  • operating_point::OperatingPoint : OperatingPoint object
  • reference_parameters::ReferenceParameters : ReferenceParameters object

Returns:

  • duct_viscous_drag::Float : total viscous drag of duct
  • boundary_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
source
DuctAPE.compute_single_side_drag_coefficientFunction
compute_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 integration
  • rotor_tip_radius::Float : radius at duct trailing edge (casing side)
  • Vref::Float : reference velocity
  • setup_boundary_layer_functions::NamedTuple : Various Akima splines and other functions for boundary layer values
  • boundary_layer_options::BoundaryLayerOptions : BoundaryLayerOptions object

Returns:

  • cd::Float : viscous drag coefficient
source
DuctAPE.compute_viscous_drag_duct_schlichtingFunction
compute_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 edge
  • duct_chord::Float : length of duct
  • TE_radius::Float : radius of the trailing edge point
  • operating_point::OperatingPoint : OperatingPoint object
source