Propulsion system

A turbofan model is provided in TASOPT.jl

Turbofan model

TASOPT.engine.tfcalc!Function
tfcalc!(wing, engine, parg, para, pare, ip, ifuel, opt_calc_call, opt_cooling, initializes_engine)

Calls on-design sizing function tfsize! or off-design analysis function tfoper! for one operating point ip.

🔃 Inputs and Outputs

Input:

  • opt_calc_call:

    • "sizing": call on-design sizing routine tfsize!
    • "oper_fixedTt4": call off-design analysis routine tfoper! with specified Tt4
    • "oper_fixedFe": call off-design analysis routine tfoper! with specified net thrust (Fe)
  • opt_cooling: turbine cooling flag

    • "none": no cooling mass flow
    • "fixed_coolingflowratio": use specified cooling flow ratios epsrow(.); calculate Tmrow(.)
    • "fixed_Tmetal": use specified metal temperatures Tmrow(.); calculate epsrow(.)
  • initializes_engine:

    • true: initialize variables for iteration in tfoper!
    • false: use current variables as initial guesses in tfoper!
source
TASOPT.engine.tfsize!Function
tfsize!(gee, M0, T0, p0, a0, M2, M25,
  Feng, Phiinl, Kinl, eng_has_BLI_cores,
  BPR, pif, pilc, pihc,
  pid, pib, pifn, pitn,
  Ttf, ifuel, etab,
  epf0, eplc0, ephc0, epht0, eplt0,
  mofft, Pofft,
  Tt9, pt9, Tt4,
  epsl, epsh,
  opt_cooling,
  Mtexit, dTstrk, StA, efilm, tfilm,
  fc0, epht_fc,
  M4a, ruc,
  ncrowx, ncrow,
  epsrow, Tmrow, 
  Δh_PreC, Δh_InterC, Δh_Regen, Δh_TurbC,
  Δp_PreC, Δp_InterC, Δp_Regen)

Turbofan performance and sizing routine.

Calculation procedure follows that of Kerrebrock, but the usual gas property formulas are replaced by function calls (described in Gas Calculations), which can therefore implement more general gas models. In addition, a turbine cooling model is added.

🔃 Inputs and Outputs

Inputs:

  • gee: gravity acceleration
  • M0: freestream Mach
  • T0: freestream temperature [K]
  • p0: freestream pressure [Pa]
  • M2: fan-face Mach number
  • M25: HPC-face Mach number
  • Feng: required net thrust (PKinl+PKout-Phi_jet)/u0 = sum( mdot u)
  • Phiinl: inlet ingested dissipation
  • eng_has_BLI_cores: false=core in clear flow, true=core sees Phiinl
  • BPR: bypass ratio = mdot_fan/mdot_core
  • pif: fan pressure ratio ( = pt7 /pt2)
  • pilc: LP comp pressure ratio ( = pt25/pt2)
  • pihc: HP comp pressure ratio ( = pt3 /pt25)
  • pid: diffuser pressure ratio ( = pt2 /pt0)
  • pib: burner pressure ratio ( = pt4 /pt3)
  • pifn: fan nozzle pressure ratio ( = pt7/pt2.1)
  • pitn: turbine nozzle pressure ratio ( = pt5/pt4.9)
  • Ttf: fuel temperature entering combustor
  • ifuel: fuel index, see function gasfun
  • hvap: fuel enthalpy of vaporization (J/kg)
  • etab: combustor efficiency (fraction of fuel burned)
  • epf0: fan max polytropic efficiency
  • eplc0: LPC max polytropic efficiency
  • ephc0: HPC max polytropic efficiency
  • epht0: HPT max polytropic efficiency
  • eplt0: LPT max polytropic efficiency
  • mofft: mass flow offtake at LPC discharge station 2.5
  • Pofft: low spool power offtake
  • Tt9: offtake air discharge total temperature
  • pt9: offtake air discharge total pressure
  • epsl: low spool power loss fraction
  • epsh: high spool power loss fraction
  • opt_cooling: turbine cooling flag
    • "none" = no cooling, ignore all cooling parameters below
    • "fixed_coolingflowratio" = usual cooling, using passed-in fcool
    • "fixed_Tmetal" = usual cooling, but set (and return) fcool from Tmetal
  • Mtexit: turbine blade-row exit Mach, for setting temperature drops
  • dTstrk: hot-streak temperature delta [K], used only if opt_cooling="fixed_Tmetal"
  • StA: area-weighted Stanton number , used only if opt_cooling="fixed_Tmetal"
  • M4a: effective Mach at cooling-flow outlet (start of mixing)
  • ruc: cooling-flow outlet velocity ratio, u/ue
  • ncrowx: dimension of epsrow array
  • ncrow: number of blade rows requiring cooling
  • epsrow(.): specified cooling-flow bypass ratio if opt_cooling="fixed_coolingflowratio"
  • Tmrow(.): specified metal temperature [K] if opt_cooling="fixed_Tmetal"

Outputs:

  • epsrow(.): resulting cooling-flow bypass ratio if opt_cooling="fixed_Tmetal"
  • Tmrow(.): resulting metal temperature [K] if opt_cooling="fixed_coolingflowratio"
  • TSFC: thrust specific fuel consumption = mdot_fuel g / F [1/s]
  • Fsp: specific thrust = F / (mdot u0) = F / ((1+BPR) mdot_core u0)
  • hfuel: fuel heating value [J / kg K]
  • ff: fuel mass flow fraction = mdot_fuel / mdot_core
  • mcore: core mass flow = mdot_core [kg/s]
  • A2: fan-face area [m^2]
  • A25: HPC-face area [m^2]
  • A5: core nozzle area [m^2]
  • A7: fan nozzle area [m^2]
  • A6: core plume area [m^2]
  • A8: fan plume area [m^2]
  • Tt?: total temperature
  • ht?: total complete enthalpy (includes heat of formation)
  • pt?: total pressure
  • cpt?: specific heat at stagnation temperature (= dh/dT)
  • Rt?: gas constant at stagnation conditions
  • T?: static temperature
  • u?: velocity
  • epf: fan polytropic efficiency
  • eplc: LPC polytropic efficiency
  • ephc: HPC polytropic efficiency
  • epht: HPT polytropic efficiency
  • eplt: LPT polytropic efficiency
  • etaf: fan overall efficiency
  • etalc: LPC overall efficiency
  • etahc: HPC overall efficiency
  • etaht: HPT overall efficiency
  • etalt: LPT overall efficiency
  • Lconv: true if convergence was successful, false otherwise

The "?" symbol denotes the station index:

  • 0: freestream
  • 18: fan face outside of casing BLs
  • 19: fan face over LPC portion
  • 2: fan face over fan portion
  • 21: fan exit
  • 25: LPC exit, HPC inlet
  • 3: compressor exit
  • 4: combustor exit before cooling air addition
  • 41: turbine inlet after cooling air addition
  • 45: HPT exit, LPT inlet
  • 49: LPT exit
  • 5: core nozzle
  • 6: core flow downstream
  • 7: fan nozzle
  • 8: fan flow downstream
source
TASOPT.engine.tfweightFunction
tfweight(ac)

Engine weight estimation function using Giulia Pantalone, Drela, or Fitzgerald model.

🔃 Inputs and Outputs

Input:

  • ac::aircraft: aircraft object

Output:

  • Weng: Total engine weight.
  • Wnac: Nacelle weight.
  • Webare: Bare engine weight.
  • Snace1: Nacelle area.
source

Turbofan Maps

TASOPT.engine.NcmapMethod
Ncmap(pratio, mb, piD, mbD, NbD, Cmap)

Calculates compressor or fan corrected speed as a function of pressure ratio and corrected mass flow

🔃 Inputs and Outputs

Inputs:

  • pratio: pressure ratio
  • mb: corrected mass flow
  • piD: design pressure ratio
  • mbD: design corrected mass flow
  • NbD: design corrected speed
  • Cmap(.): map constants

Outputs:

  • Nb: wheel speed
  • Nb_?: derivatives
source
TASOPT.engine.ecmapMethod
ecmap(pratio, mb, piD, mbD, Cmap, effo, piK, effK)

Calculates compressor or fan efficiency as a function of pressure ratio and corrected mass flow

🔃 Inputs and Outputs

Inputs:

  • pratio: pressure ratio
  • mb: corrected mass flow
  • piD: design pressure ratio
  • mbD: design corrected mass flow
  • Cmap(.): map constants
  • effo: maximum efficiency
  • piK: pi-dependence offset eff = effo + effK*(pi-piK)
  • effK: pi-dependence slope

Outputs:

  • eff: efficiency
  • eff_?: derivatives
source
TASOPT.engine.Ncmap1Method
Ncmap1(pratio, m, piD, mbD, NbD, ABCDm, iabcd, Tr, pr)

Calculates compressor or fan efficiency as a function of pressure ratio and corrected mass flow

🔃 Inputs and Outputs

Inputs:

  • pratio: pressure ratio
  • mb: corrected mass flow
  • piD: design pressure ratio
  • mbD: design corrected mass flow
  • NbD: design corrected speed
  • ABCDm: map constants
  • iabcd: map exponents
  • Tr: T/Tref
  • pr: p/pref

Outputs:

  • N: wheel speed
  • N_?: derivatives
source
TASOPT.engine.ecmap1Method
ecmap1(pratio, m, piD, mbD, ABCDm, iabcd, effo, Tr, pr)

Calculates compressor or fan efficiency as a function of pressure ratio and mass flow

🔃 Inputs and Outputs

Inputs:

  • pratio: pressure ratio
  • mb: corrected mass flow
  • piD: design pressure ratio
  • mbD: design corrected mass flow
  • ABCDm: map constants
  • iabcd: map exponents
  • effo: maximum efficiency
  • Tr: T/Tref
  • pr: p/pref

Outputs:

  • eff: efficiency
  • eff_?: derivatives
source
TASOPT.engine.etmapMethod
etmap(dh, mb, Nb, piD, mbD, NbD, ept0, Tmap, Tt, cpt, Rt)

Calculates turbine efficiency as a function of work and corrected mass flow

🔃 Inputs and Outputs

Inputs:

  • dh: enthalpy change
  • mb: corrected mass flow
  • Nb: corrected speed
  • piD: design pressure ratio
  • mbD: design corrected mass flow
  • NbD: design corrected speed
  • ept0: turbine polytropic efficiency estimate
  • Tmap(.): map constants

Outputs:

  • eff: efficiency
  • eff_?: derivatives
source
TASOPT.engine.PimapMethod
Pimap(mb, Nb, piD, mbD, NbD, Cmap)

Calculates compressor or fan pressure ratio as a function of pressure ratio and corrected mass flow

🔃 Inputs and Outputs

Inputs:

  • mb: corrected mass flow
  • Nb: corrected speed
  • piD: design pressure ratio
  • mbD: design corrected mass flow
  • NbD: design corrected speed
  • Cmap(.): map constants

Outputs:

  • pratio: pressure ratio
  • pi_?: derivatives
source
TASOPT.engine.tfoper!Function
function tfoper!(gee, M0, T0, p0, a0, Tref, pref,
  Phiinl, Kinl, eng_has_BLI_cores,
  pid, pib, pifn, pitn,
  Gearf,
  pifD, pilcD, pihcD, pihtD, piltD,
  mbfD, mblcD, mbhcD, mbhtD, mbltD,
  NbfD, NblcD, NbhcD, NbhtD, NbltD,
  A2, A25, A5, A7,
  opt_calc_call,
  Ttf, ifuel, etab,
  epf0, eplc0, ephc0, epht0, eplt0,
  mofft, Pofft,
  Tt9, pt9,
  epsl, epsh,
  opt_cooling,
  Mtexit, dTstrk, StA, efilm, tfilm,
  M4a, ruc,
  ncrowx, ncrow,
  epsrow, Tmrow,
  Feng,
  M2, pif, pilc, pihc, mbf, mblc, mbhc, Tt4, pt5, mcore, M25)

Turbofan operation routine. Calculation procedure follows that of Kerrebrock, but the usual gas property formulas are replaced by function calls, which can therefore implement more general gas models. In addition, a turbine cooling model is added.

The gas routines are described in Gas Calculations

🔃 Inputs and Outputs

Inputs:

  • gee: gravity acceleration

  • M0: freestream Mach

  • T0: freestream temperature [K]

  • p0: freestream pressure [Pa]

  • Tref: reference temperature for corrected mass flow and speed

  • pref: reference pressure for corrected mass flow

  • Phiinl: inlet ingested dissipation Phi_inl

  • eng_has_BLI_cores:

    • false: core in clear flow
    • true: core sees Phiinl
  • pid: diffuser pressure ratio ( = pt2/pt0)

  • pib: burner pressure ratio ( = pt4/pt3)

  • pifn: fan nozzle pressure ratio ( = pt7/pt6.9)

  • pitn: turbine nozzle pressure ratio ( = pt5/pt4.9)

  • Gearf: fan gear ratio ( = Nl/Nf )

  • pifD: design fan pressure ratio ( = pt21/pt2)

  • pilcD: design LPC pressure ratio ( = pt25/pt19)

  • pihcD: design HPC pressure ratio ( = pt3/pt25)

  • pihtD: design HPT pressure ratio ( = pt45/pt41)

  • piltD: design LPT pressure ratio ( = pt49/pt45)

  • mbfD: design corrected fan mass flow ( = mf*sqrt(Tt2/Tref)/(pt2/pref) ) where mf = mc*BPR

  • mblcD: design corrected LPC mass flow ( = mc*sqrt(Tt19/Tref)/(pt19/pref) )

  • mbhcD: design corrected HLC mass flow ( = mc*sqrt(Tt25/Tref)/(pt25/pref) )

  • mbhtD: design corrected HPT mass flow ( = mt*sqrt(Tt41/Tref)/(pt41/pref) ) where mt = mc*(1+ff)

  • mbltD: design corrected LPT mass flow ( = mt*sqrt(Tt45/Tref)/(pt45/pref) )

  • NbfD: design corrected fan speed ( = Nf/sqrt(Tt2/Tref) )

  • NblcD: design corrected LPC speed ( = Nl/sqrt(Tt19/Tref) )

  • NbhcD: design corrected HPC speed ( = Nh/sqrt(Tt25/Tref) )

  • NbhtD: design corrected HPT speed ( = Nh/sqrt(Tt41/Tref) )

  • NbltD: design corrected LPT speed ( = Nl/sqrt(Tt45/Tref) )

  • A2: fan-face area [m^2]

  • A25: HPC-face area [m^2]

  • A5: core nozzle area [m^2]

  • A7: fan nozzle area [m^2]

  • opt_calc_call:

    • "oper_fixedTt4": Tt4 is specified
    • "oper_fixedFe": Feng is specified
  • Tt4: turbine-inlet total temperature [K]

  • Ttf: fuel temperature entering combustor

  • ifuel: fuel index, see function gasfun

  • hvap: fuel enthalpy of vaporization (J/kg)

  • etab: combustor efficiency (fraction of fuel burned)

  • epf0: max fan polytropic efficiency

  • eplc0: LPC max polytropic efficiency

  • ephc0: HPC max polytropic efficiency

  • epht0: HPT max polytropic efficiency

  • eplt0: LPT max polytropic efficiency

  • mofft: mass flow offtake at LPC discharge station 2.5

  • Pofft: low spool power offtake

  • Tt9: offtake air discharge total temperature

  • pt9: offtake air discharge total pressure

  • epsl: low spool power loss fraction

  • epsh: high spool power loss fraction

  • opt_cooling: turbine cooling flag

    • "none": no cooling, ignore all cooling parameters below
    • "fixed_coolingflowratio": usual cooling, using passed-in fcool
    • "fixed_Tmetal": usual cooling, but set (and return) fcool from Tmetal
  • Mtexit: turbine blade-row exit Mach, for setting temperature drops

  • Tmetal: specified metal temperature [K], used only if opt_cooling="fixed_Tmetal"

  • dTstrk: hot-streak temperature delta {K}, used only if opt_cooling="fixed_Tmetal"

  • StA: area-weighted Stanton number , used only if opt_cooling="fixed_Tmetal"

  • M4a: effective Mach at cooling-flow outlet (start of mixing)

  • ruc: cooling-flow outlet velocity ratio, u/ue

  • ncrowx: dimension of epsrow array

  • ncrow: number of blade rows requiring cooling

  • epsrow(.): input specified cooling-flow bypass ratio if opt_cooling="fixed_coolingflowratio"; output resulting cooling-flow bypass ratio if opt_cooling="fixed_Tmetal".

  • Tmrow(.): input specified metal temperature [K] if opt_cooling="fixed_Tmetal"; output resulting metal temperature [K] if opt_cooling="fixed_coolingflowratio"

Output:

  • epsrow(.): see above
  • Tmrow(.): see above
  • TSFC: thrust specific fuel consumption = mdot_fuel g / F [1/s]
  • Fsp: specific thrust = F / (mdot u0) = F / ((1+BPR) mdot_core u0)
  • hfuel: fuel heating value [J / kg K]
  • ff: fuel mass flow fraction = mdot_fuel / mdot_core
  • Feng: net effective thrust = (PK_inl+PK_out-Phi_jet)/u0 = sum(mdot u)
  • mcore: core mass flow = mdot_core [kg/s]
  • BPR: bypass ratio = mdot_fan/mdot_core
  • Tt?: total temperature
  • ht?: total complete enthalpy (includes heat of formation)
  • pt?: total pressure
  • cpt?: specific heat at stagnation temperature (= dh/dT)
  • Rt?: gas constant at stagnation conditions
  • T?: static temperature
  • u?: velocity
  • etaf: fan overall efficiency
  • etac: compressor overall efficiency
  • etatf: fan-turbine overall efficiency
  • etatc: comp-turbine overall efficiency
  • Lconv: true if convergence was successful, false otherwise

The "?" symbol denotes the station index:

  • 0: freestream
  • 18: fan face outside of casing BLs
  • 19: fan face over LPC portion
  • 2: fan face over fan portion
  • 21: fan exit, precooler inlet
  • 19c: precooler outlet, LPC inlet
  • 25: LPC exit, intercooler inlet
  • 25c: intercooler exit, HPC inlet
  • 3: compressor exit
  • 4: combustor exit before cooling air addition
  • 41: turbine inlet after cooling air addition
  • 45: HPT exit, LPT inlet
  • 49: LPT exit, regenerative cooler inlet
  • 49c: regenerative cooler outlet
  • 5: core nozzle
  • 6: core flow downstream
  • 7: fan nozzle
  • 8: fan flow downstream
source

Turbofan Cooling

TASOPT.engine.mcoolMethod
mcool(ncrowx, Tmrow, Tt3, Tt4, dTstreak, Trrat, efilm, tfilm, StA)

Calculates cooling mass flow requirement.

🔃 Inputs and Outputs

Inputs:

  • ncrowx: dimension of Tmrow(.) and epsrow(.) arrays (max number of blade rows)
  • Tmrow(.): design metal temperature for each blade row
  • Tt3: cooling flow temperature
  • Tt4: hot gas temperature from burner
  • dTstreak: hot-streak temperature increase over Tt4, for first blade row
  • Trrat: static temperature ratio across each blade row, T4.1 / T4
  • efilm: cooling efficiency = (Tco-Tci)/(Tmetal-Tci)
  • tfilm: film effectiveness = (Tgas-Tfaw)/(Tgas-Tco) where
    • Tco = temperature of cooling air exiting blade
    • Tci = temperature of cooling air entering blade
    • Tfaw = film adiabatic wall temperature (for insulated-wall case)
  • StA: area-weighted external Stanton number =St (Asurf/Aflow) cpgas/cpcool`

Output:

  • ncrow: number of blade rows which need cooling
  • epsrow(.): cooling mass flow ratio for each blade row, m_c_row/m_air
source
TASOPT.engine.TmcalcMethod
Tmcalc(ncrowx, ncrow, Tt3, Tt4, dTstreak, Trrat, efilm, tfilm, StA, epsrow)

Calculates metal temperature for blade row

🔃 Inputs and Outputs

Inputs:

  • ncrowx: dimension of Tmrow(.),epsrow(.) arrays (max number of blade rows)
  • ncrow: number of blade rows which are cooled
  • epsrow(.): cooling mass flow ratio for each blade row, mcrow/m_air
  • Tt3: cooling flow temperature
  • Tt4: hot gas temperature from burner
  • dTstreak: hot-streak temperature increase over Tt4, for first blade row
  • Trrat: static temperature ratio across each blade row, T4.1 / T4
  • efilm: cooling efficiency = (Tco-Tci)/(Tmetal-Tci)
  • tfilm: film effectiveness = (Tgas-Tfaw)/(Tgas-Tco) Tco = temperature of cooling air exiting blade Tci = temperature of cooling air entering blade Tfaw = film adiabatic wall temperature (for insulated-wall case)
  • StA: area-weighted external Stanton number = St (Asurf/Aflow) cpgas/cpcoolOutput:
  • Tmrow(.): design metal temperature for each blade row
source

Turbomachinery Components

The compressor off-design performance is determined by interpolation to the compressor maps in pyCycle[1]. The compressor parameters are scaled to the design pressure ratios, speeds, and mass flow rates in the pyCycle maps by using

\[ \tilde{p} = \frac{\pi -1}{\pi_D -1}\]

\[ \tilde{m} = \frac{\bar{m}}{\bar{m}_D}\]

\[ \tilde{N} = \frac{\bar{N}}{\bar{N}_D},\]

where $\pi$ represents the pressure ratio, $\bar{m}$ is the corrected mass flow rate, $\bar{N}$ is the corrected speed, and the subscript $D$ denotes the design values. For a given set of compressor parameters, the normalized parameters $\tilde{p}$, $\tilde{m}$, and $\tilde{N}$ are used to calculate the dimensional values in the pyCycle map space.

The pyCycle maps contain data for pressure ratio, corrected mass flow rate, and isentropic efficiency as a function of corrected speed and R-line parameter. However, the turbofan operation function, tfoper!(), is set up so that the polytropic efficiency and corrected speed is calculated from the pressure ratio and corrected mass flow rate. Therefore, a reverse interpolation problem is required to compute these parameters. For this, the compressor maps are extrapolated and a standard non-linear solver is used in find_NR_inverse_with_derivatives() to calculate the corrected speed and R-line parameter that correspond to given corrected mass flow rate and pressure ratio. The extrapolated pressure ratio and polytropic efficiency maps are shown below. These parameters are then translated to the scaled map by using $\tilde{p}$, $\tilde{m}$, and $\tilde{N}$ in calculate_compressor_speed_and_efficiency(), which also returns the polytropic efficiency and derivatives.

PEMfig PEMfig

TASOPT.engine.calculate_compressor_speed_and_efficiencyFunction
calculate_compressor_speed_and_efficiency(map, pratio, mb, piD, mbD, NbD, epol0; Ng=0.5, Rg=2.0)

Calculates corrected speed and polytropic efficiency for a compressor, along with derivatives to pressure ratio and mass flow.

🔃 Inputs and Outputs

Inputs:

  • map::CompressorMap: structure containing the compressor map and interpolations.
  • pratio::Float64: compressor pressure ratio.
  • mb::Float64: current mass flow rate.
  • piD::Float64: design pressure ratio.
  • mbD::Float64: design mass flow rate.
  • NbD::Float64: design rotational speed.
  • epol0::Float64: maximum polytropic efficiency.
  • Ng::Float64: initial guess for normalized speed (optional).
  • Rg::Float64: initial guess for R-line (optional).

Outputs:

  • Nb::Float64: corrected compressor speed.
  • epol::Float64: polytropic efficiency.
  • dNb_dpi::Float64: derivative of corrected speed w.r.t. piD.
  • dNb_dmb::Float64: derivative of corrected speed w.r.t. mb.
  • depol_dpi::Float64: derivative of efficiency w.r.t. piD.
  • depol_dmb::Float64: derivative of efficiency w.r.t. mb.
  • N::Float64: matched normalized speed.
  • R::Float64: matched R-line.
source
TASOPT.engine.find_NR_inverse_with_derivativesFunction
find_NR_inverse_with_derivatives(itp_Wc, itp_PR, Wc_target, PR_target; Ng=0.5, Rg=2.0)

Finds the normalized speed N and R-line R corresponding to a target corrected mass flow Wc_target and pressure ratio PR_target, using a nonlinear solver with Jacobian information from interpolation gradients.

🔃 Inputs and Outputs

Inputs:

  • itp_Wc::GriddedInterpolation: interpolation of corrected mass flow over (N, R).
  • itp_PR::GriddedInterpolation: interpolation of pressure ratio over (N, R).
  • Wc_target::Float64: target corrected mass flow.
  • PR_target::Float64: target pressure ratio.
  • Ng::Float64: initial guess for normalized speed (optional).
  • Rg::Float64: initial guess for R-line (optional).

Outputs:

  • N::Float64: normalized speed at the matched point.
  • R::Float64: R-line value at the matched point.
  • dN_dw::Float64: ∂N/∂Wc
  • dN_dpr::Float64: ∂N/∂PR
  • dR_dw::Float64: ∂R/∂Wc
  • dR_dpr::Float64: ∂R/∂PR
source
  • 1https://github.com/OpenMDAO/pyCycle