Welcome to Bamboo’s documentation!

Use the sidebar (on the left) to find lists of all the available classes and functions, for each submodule. Use the ‘home’ buttons (with the small image of a house) to navigate back to this main page.

cusfbamboo.engine

Main Engine class, as well as tools for specifying engine geometry and the perfect gas model used to calculate flow properties.

cusfbamboo.plot

Module that provides plotting tools, to streamline the creation of plots.

cusfbamboo.materials

Classes for representing solid materials, fluids and transport properties.

cusfbamboo.rao

Rao bell nozzle data.

cusfbamboo.isen

Isentropic compressible flow relations.

cusfbamboo.hx

General solver for coflow and counterflow heat exchangers, using a 1-D thermal resistance model.

cusfbamboo.circuit

Classes and functions related to thermal circuit calculations.

cusfbamboo.engine module

Main Engine class, as well as tools for specifying engine geometry and the perfect gas model used to calculate flow properties.

References

Notes

  • With some exceptions, coolant properties are currently evaluated at the bulk temperature, instead of the film temperature. This is because the high wall temperature can sometimes be above the coolant boiling point, which can cause errors. Ideally you would use nucleate boiling correlations in this case, but this is not always possible, and so use the bulk properties is available as a compromise.

class cusfbamboo.engine.PerfectGas(**kwargs)

Bases: object

Object to store a perfect gas model (i.e. an ideal gas with constant cp and cv). You only need to input 2 properties to fully define it.

Keyword Arguments
  • gamma (float) – Ratio of specific heats cp/cv.

  • cp (float) – Specific heat capacity at constant pressure (J/kg/K)

  • molecular_weight (float) – Molecular weight of the gas (kg/kmol)

gamma

Ratio of specific heats cp/cv.

Type

float

cp

Specific heat capacity at constant pressure (J/kg/K)

Type

float

molecular_weight

Molecular weight of the gas (kg/kmol)

Type

float

R

Specific gas constant (J/kg/K)

Type

float

class cusfbamboo.engine.ChamberConditions(p0, T0)

Bases: object

Object for storing combustion chamber thermodynamic conditions. The mass flow rate does not have to be defined - it is fixed by the nozzle throat area.

Parameters
  • p0 (float) – Gas stagnation pressure (Pa).

  • T0 (float) – Gas stagnation temperature (K).

class cusfbamboo.engine.Geometry(xs, rs)

Bases: object

Class for representing the inner contour of a rocket engine, from the beginning of the combustion chamber to the nozzle exit.

Parameters
  • xs (list) – Array of x-positions, that the ‘y’ list corresponds to (m). Must be increasing values of x.

  • rs (list) – Array, containing local engine radius (m).

x_t

x-position of the throat (m)

Type

float

r_t

Throat radius (m)

Type

float

A_t

Throat area (m2)

Type

float

r_e

Exit radius (m)

Type

float

A_e

Exit area (m2)

Type

float

r_curvature_t

Radius of curvature at the throat (m)

Type

float

property x_t
property r_t
property A_t
property r_e
property A_e
plot()

Plot the engine geometry. Must run cusfbamboo.plot.show() or matplotlib.pyplot.show() to see the plot.

r(x)

Get the distance from the centreline to the inner wall of the engine.

Parameters

x (float) – x position (m)

Returns

Distance from engine centreline to edge of inner wall (m)

Return type

float

dr_dx(x)

Get the slope of the engine wall, dr/dx.

Parameters

x (float) – Axial position (m).

Returns

Rate of change of contour radius with respect to position, dr/dx

Return type

float

A(x)

Get the flow area for the exhaust gas

Parameters

x (float) – x position (m)

Returns

Flow area (m2)

Return type

float

class cusfbamboo.engine.Wall(material, thickness)

Bases: object

Object for representing an engine wall.

Parameters
  • material (Material) – Material object to define the material the wall is made of.

  • thickness (float or callable) – Thickness of the wall (m). Can be a constant float, or a function of position, i.e. t(x).

thickness(x)

Get the thickness of the wall at a position x.

Parameters

x (float) – Axial position along the engine (m)

Returns

Wall thickness (m)

Return type

float

class cusfbamboo.engine.CoolingJacket(T_coolant_in, p_coolant_in, mdot_coolant, channel_height, coolant_transport, roughness=None, configuration='vertical', **kwargs)

Bases: object

Class for representing cooling jacket properties.

Note

Spiralling channels are assumed to cover the entire surface area of the outer chamber wall. A blockage ratio can still be used to ‘block up’ part of the channels with fins.

Note

All channels are assumed to have a rectangular cross section.

Parameters
  • T_coolant_in (float) – Inlet static temperature of the coolant (K)

  • p_coolant_in (float) – Inlet static pressure of the coolant (Pa)

  • mdot_coolant (float) – Mass flow rate of the coolant (kg/s)

  • channel_height (float or callable) – Radial height of the cooling channels (i.e. the distance between the inner and outer wall that the coolant flows through) (m). Can be a constant float, or function of axial position (x).

  • coolant_transport (TransportProperties) – Transport properties of the coolant

  • configuration (str, optional) – Type of cooling channel. Either ‘vertical’ for straight, axial, channels. Or ‘spiral’ for a helix around the engine. Defaults to “vertical”.

  • roughness (float or callable, optional) – Wall roughness in the channel (m), for pressure drop calculations. Can be a constant float, or function of axial position (x). Defaults to None, in which case a smooth walled approximation is used.

Keyword Arguments
  • blockage_ratio (float or callable) – This is the proportion (by area) of the channel cross section occupied by fins. Can be a constant float, or a function of axial position (x). Defaults to zero.

  • number_of_channels (int) – The number of parallel cooling channels (with fins in between them). Only relevant if ‘blockage_ratio’ !=0.

  • channel_width (float or callable) – Width of a single coolant channel. Can be a constant float, or a function of axial position (x).

  • xs (list) – Minimum and maximum x value which the cooling jacket is present over (m), e.g. (x_min, x_max). Can be in either order.

  • restrain_fins (bool) – Whether or not the fins in cooling channels are physically restrained by (i.e. attached to) the outer cooling jacket. This affects the pressure stress. Automatically ignored if blockage_ratio = 0. Defaults to False.

channel_height(x)

Get the channel height at a position, x.

Parameters

x (float) – Axial position along the engine (m)

Returns

Channel height (m)

Return type

float

blockage_ratio(x)

Get the blockage ratio at a position, x.

Parameters

x (float) – Axial position along the engine (m)

Returns

Blockage ratio

Return type

float

channel_width(x)

Get the width of a single cooling channel for a spiral configuration, at a position, x.

Parameters

x (float) – Axial position along the engine (m)

Returns

Channel width (m)

Return type

float

bundle_width(x)

Width of a ‘bundle’ of spiralling cooling channels. Equal to channel_width * number_of_channels.

Parameters

x (float) – Axial position along the engine (m)

Returns

Width per channel multiplied by number of channels (m)

Return type

float

roughness(x)

Get the channel roughness, at a position, x.

Parameters

x (float) – Axial position along the engine (m)

Returns

Wall roughness of the channel (m)

Return type

float

f_darcy_laminar(ReDh, Dh, x)
f_darcy_turbulent(ReDh, Dh, x)
f_darcy(ReDh, Dh, x)
class cusfbamboo.engine.Engine(perfect_gas, chamber_conditions, geometry, coolant_convection='gnielinski', exhaust_convection='bartz-sigma', **kwargs)

Bases: object

Class for representing a liquid rocket engine.

Parameters
  • perfect_gas (PerfectGas) – PerfectGas representing the exhaust gas for the engine.

  • chamber_conditions (ChamberConditions) – ChamberConditions for the engine.

  • geometry (Geometry) – Geomtry object to define the engine’s contour.

  • coolant_convection (str) – Convective heat transfer model to use for the coolant side. Can be ‘dittus-boelter’, ‘sieder-tate’ or ‘gnielinski’. Defaults to ‘gnielinski’.

  • exhaust_convection (str) – Convective heat transfer model to use the for exhaust side. Can be ‘dittus-boelter’, ‘bartz’ or ‘bartz-sigma’. Defaults to ‘bartz-sigma’.

Keyword Arguments
  • walls (Wall or list) – Either a single Wall object that specifies the combustion chamber wall, or a list of Wall objects that represent multiple layers with different materials. List must be in the order [hottest_wall, … , coldest_wall].

  • cooling_jacket (CoolingJacket) – CoolingJacket object to specify the cooling jacket on the engine.

  • exhaust_transport (TransportProperties) – TransportProperties object that defines the exhaust gas transport properties.

mdot

Mass flow rate of exhaust gas (kg/s)

Type

float

c_star

C* for the engine (m/s).

Type

float

coolant_convection

Convective heat transfer model to use for the coolant side.

Type

str

exhaust_convection

Convective heat transfer model to use the for exhaust side.

Type

str

h_exhaust_sf

Scale factor for the exhaust convective heat transfer coeffcient. Defaults to 1.

Type

float

h_coolant_sf

Scale factor for the coolant convective heat transfer coeffcient. Defaults to 1.

Type

float

walls

List of Wall objects between the hot gas and coolant

Type

list

M(x)

Get exhaust gas Mach number.

Parameters

x (float) – Axial position along the engine (m).

Returns

Mach number of the freestream.

Return type

float

T(x)

Get temperature at a position along the nozzle. :param x: Distance from the throat, along the centreline (m) :type x: float

Returns

Temperature (K)

Return type

float

p(x)

Get pressure at a position along the nozzle. :param x: Distance from the throat, along the centreline (m) :type x: float

Returns

Pressure (Pa)

Return type

float

rho(x)

Get exhaust gas density. :param x: Axial position. Throat is at x = 0. :type x: float

Returns

Freestream gas density (kg/m3)

Return type

float

total_wall_thickness(x)
plot()

Plot the engine geometry, including the cooling channels, all to scale. You will need to run matplotlib.pyplot.show() or cusfbamboo.plot.show() to see the plot.

helix_angle(x)

Angle between the spiralling helix and the axial direction. Same as beta in Ref [7].

Parameters

x (float) – Axial position (m)

Returns

Helix angle (rad)

Return type

float

coolant_pitch(x)
coolant_slope(x)

Angle that the coolant wall is sloped at radially, same as alpha in Ref [7]

Parameters

x (float) – Axial position (m)

Returns

Coolant wall slope (rad)

Return type

float

dLc_dx(x)

Conversion factor between travelling a distance ‘dx’ in the axial direction, and the path taken by the cooling fluid ‘dLc’.

Parameters

x (float) – Axial position x (m).

Returns

dLc/dx

Return type

float

A_coolant(x)

Flow area of the coolant at an axial position.

Parameters

x (float) – Axial position x (m)

Returns

Coolant flow area (m2)

Return type

float

Dh_coolant(x)

Hydraulic diameter of the coolant flow channel - used for pressure drops. This is equal to 4 * A / P, where ‘A’ is the coolant flow area and ‘P’ is the perimeter of the channel.

Parameters

x (float) – Axial position (m)

Returns

Hydraulic diameter (m)

Return type

float

thrust(p_amb=100000.0)

Get the thrust of the engine for a given ambient pressure

Parameters

p_amb (float, optional) – Ambient pressure (Pa). Defaults to 1e5.

Returns

Thrust (N)

Return type

float

isp(p_amb=100000.0)

Get the specific impulse of the engine for a given ambient pressure.

Parameters

p_amb (float, optional) – Ambient pressure (Pa). Defaults to 1e5.

Returns

Specific impulse (m/s)

Return type

float

T_h(state)
cp_c(state)
A_c(state)
V_c(state)
Rdx(state)
extra_dQ_dx(state)
dp_dx_f(state)
steady_heating_analysis(num_grid=1000, counterflow=True, iter_start=5, iter_each=2)

Run a steady state cooling simulation.

Parameters
  • num_grid (int) – Number of grid points to use (1-dimensional)

  • counterflow (bool, optional) – Whether or not the cooling is flowing coutnerflow or coflow, relative to the exhaust gas. Defaults to True (which means counterflow).

  • iter_start (int) – Number of times to iterate on the entry conditions. Defaults to 5.

  • iter_each (int) – Number of times to iterate on the solution at each datapoint. Defaults to 2.

cusfbamboo.plot module

Module that provides plotting tools, to streamline the creation of plots.

cusfbamboo.plot.show()

Show a plot. Internally, this simply runs matplotlib.pyplot.show().

cusfbamboo.plot.plot_temperatures(data_dict, only_indexes=None)

Given the output dictionary from an engine cooling analysis, plot the temperatures against position. Note you will have to run matplotlib.pyplot.show() or cusfbamboo.plot.show() to see the plot.

Parameters
  • data_dict (dict) – Dictionary contaning the cooling analysis results.

  • only_indexes (list) – List of temperature indexes to plot, e.g. [1, -1] would only show the coolant temperature and exhaust temperature. Defaults to None, which means everything is plotted.

cusfbamboo.plot.plot_p_coolant(data_dict)

Given the output dictionary from a engine cooling analysis, plot the cooling jacket static pressure against x position. Note you will have to run matplotlib.pyplot.show() or cusfbamboo.plot.show() to see the plot.

Parameters

data_dict (dict) – Dictionary contaning the cooling analysis results.

cusfbamboo.plot.plot_T_coolant(data_dict)

Given the output dictionary from a engine cooling analysis, plot the coolant static temperature against postiion

Parameters

data_dict (dict) – Dictionary contaning the cooling analysis results.

cusfbamboo.plot.plot_q_per_area(data_dict)

Given the output dictionary from a engine cooling analysis, plot the heat flux against position. Note you will have to run matplotlib.pyplot.show() or cusfbamboo.show() to see the plot.

Parameters

data_dict (dict) – Dictionary contaning the cooling analysis results.

cusfbamboo.plot.plot_tangential_stress(data_dict, wall_index=0)

Given the output dictionary from a engine cooling analysis, plot the thermal stress in the inner chamber wall, against position. Note you will have to run matplotlib.pyplot.show() or bam.show() to see the plot.

Parameters
  • data_dict (dict) – Dictionary contaning the cooling analysis results.

  • wall_index (int) – The index of the wall to plot the stresses for. Defaults to 0 (the exhaust side wall).

cusfbamboo.plot.plot_coolant_velocity(data_dict)

Given the output dictionary from a engine cooling analysis, plot the cooling velocity against axial position. Note you will have to run matplotlib.pyplot.show() or cusfbamboo.show() to see the plot.

Parameters

data_dict (dict) – Dictionary contaning the cooling analysis results.

cusfbamboo.plot.plot_coolant_density(data_dict)

Plot the cooling density against axial position. Note you will have to run matplotlib.pyplot.show() or cusfbamboo.show() to see the plot.

Parameters

data_dict (dict) – Dictionary contaning the cooling analysis results.

cusfbamboo.plot.plot_thermal_resistances(data_dict, only_indexes=None)

Given the output dictionary from an engine cooling analysis, plot the local thermal resistances against position. Note you will have to run matplotlib.pyplot.show() or cusfbamboo.plot.show() to see the plot.

Parameters
  • data_dict (dict) – Dictionary contaning the cooling analysis results.

  • only_indexes (list) – List of resistance indexes to plot, e.g. [1, -1] would only show the coolant convective resistance and exhaust convective resistance. Defaults to None, which means everything is plotted.

cusfbamboo.plot.plot_coolant_h(data_dict)

Plot convective heat transfer coefficient for the coolant side.

Parameters

data_dict (dict) – Dictionary contaning the cooling analysis results.

cusfbamboo.materials module

Classes for representing solid materials, fluids and transport properties. Also contains some useful pre-defined materials.

References *need to get solid materials references - [1] - CoolProp, http://coolprop.org/

class cusfbamboo.materials.Material(k, **kwargs)

Bases: object

Class used to specify a material and its properties. For calculating temperatures, only ‘k’ must be defined. For stresses, you also need E, alpha, and poisson.

Parameters

k (float) – Thermal conductivity (W/m/K)

Keyword Arguments
  • E (float) – Young’s modulus (Pa)

  • alpha (float) – Thermal expansion coefficient (strain/K)

  • poisson (float) – Poisson’s ratio

class cusfbamboo.materials.TransportProperties(Pr, mu, k, cp=None, rho=None, gamma_coolant=None)

Bases: object

Container for specifying your transport properties. Each input can either be a function of temperature (K) and pressure (Pa) in that order, e.g. mu(T, p). Otherwise they can be constant floats.

Parameters
  • Pr (float or callable) – Prandtl number.

  • mu (float or callable) – Absolute viscosity (Pa s).

  • k (float or callable) – Thermal conductivity (W/m/K).

  • cp (float or callable, optional) – Isobaric specific heat capacity (J/kg/K) - only required for coolants.

  • rho (float or callable, optional) – Density (kg/m^3) - only required for coolants.

  • gamma_coolant (float or callable, optional) – Ratio of specific heats (cp/cv) for a compressible coolant. If this is submitted, it is assumed that this object represents a compressible coolant.

compressible_coolant

Whether or not this TransportProperties object represents a compressible coolant.

Type

bool

Pr(T, p)

Prandtl number.

Parameters
  • T (float) – Temperature (K)

  • p (float) – Pressure (Pa)

Returns

Prandtl number

Return type

float

mu(T, p)

Absolute viscosity (Pa s)

Parameters
  • T (float) – Temperature (K)

  • p (float) – Pressure (Pa)

Returns

Absolute viscosity (Pa s)

Return type

float

k(T, p)

Thermal conductivity (W/m/K)

Parameters
  • T (float) – Temperature (K)

  • p (float) – Pressure (Pa)

Returns

Thermal conductivity (W/m/K)

Return type

float

rho(T, p)

Density (kg/m^3) :param T: Temperature (K) :type T: float :param p: Pressure (Pa) :type p: float

Returns

Density (kg/m^3)

Return type

float

cp(T, p)

Isobaric specific heat capacity (J/kg/K)

Parameters
  • T (float) – Temperature (K)

  • p (float) – Pressure (Pa)

Returns

Isobaric specific heat capacity (J/kg/K)

Return type

float

gamma_coolant(T, p)

Ratio of specific heat capacities for a compressible coolant.

Parameters
  • T (float) – Temperature (K)

  • p (float) – Pressure (Pa)

Returns

Ratio of specific heat capacities (cp/cv).

Return type

float

class cusfbamboo.materials.NucleateBoiling(vapour_transport, liquid_transport, sigma, h_fg, C_sf)

Bases: object

Class for representing the information needed to model nucleate boiling. Not currently used.

Parameters
  • vapour_transport (TransportProperties) – The transport properties of the vapour phase.

  • liquid_transport (TransportProperties) – The transport properties of the liquid phase.

  • sigma (callable) – Surface tension of the liquid-vapour interface (N/m), as a function of temperature (K) and pressure (Pa) in the form sigma(T,p).

  • h_fg (callable) – Enthalpy between vapour and liquid phases, as a function of pressure (Pa). h_fg = h_g - h_f. (J/kg/K)

  • C_sf (float) – Surface-fluid coefficient. Will be different for different material + fluid combinations. Some examples are available in References [4] and [6] given in bamboo.circuit.py.

cusfbamboo.rao module

Rao bell nozzle data.

References

cusfbamboo.rao.rao_theta_n(area_ratio, length_fraction=0.8)

Returns the contour angle at the inflection point of the bell nozzle, by interpolating data. Data obtained by using http://www.graphreader.com/ on the graph in Reference [1].

Parameters
  • area_ratio (float) – Area ratio of the nozzle (A2/At)

  • length_fraction (int, optional) – Nozzle contraction percentage, as defined in Reference [1]. Defaults to 0.8.

Returns

“theta_n”, angle at the inflection point of the bell nozzle (rad)

Return type

float

cusfbamboo.rao.rao_theta_e(area_ratio, length_fraction=0.8)

Returns the contour angle at the exit of the bell nozzle, by interpolating data. Data obtained by using http://www.graphreader.com/ on the graph in Reference [1].

Parameters
  • area_ratio (float) – Area ratio of the nozzle (A2/At)

  • length_fraction (int, optional) – Nozzle contraction percentage, as defined in Reference [1]. Defaults to 0.8.

Returns

“theta_e”, angle at the exit of the bell nozzle (rad)

Return type

float

cusfbamboo.rao.get_rao_contour(r_c, r_t, area_ratio, L_c, theta_conv=45)

Get the x and y positions for an 80% Rao bell nozzle

Parameters
  • r_c (float) – Chamber radius (m)

  • r_t (float) – Throat radius (m)

  • area_ratio (float) – The area ratio (exit area / throat area).

  • L_c (float) – Chamber length, from the injector to the start of the nozzle converging section (m)

  • theta_conv (int, optional) – Angle of converging nozzle section (deg). Defaults to 45.

Returns

Nozzle x coordinates and corresponding radii, xs, rs (m)

Return type

(list, list)

cusfbamboo.isen module

Isentropic compressible flow relations.
cusfbamboo.isen.m_bar(M, gamma)

Non-dimensional mass flow rate, defined as m_bar = mdot * sqrt(cp*T0)/(A*p0). A is the local cross sectional area that the flow is moving through.

Parameters
  • M (float) – Mach number

  • gamma (float) – Ratio of specific heats cp/cv

Returns

Non-dimensional mass flow rate

Return type

float

cusfbamboo.isen.A_At(M, gamma)

Ratio of local flow area to the throat area, for a given Mach number [1].

Parameters
  • M (float) – Mach number

  • gamma (float) – Ratio of specific heats cp/cv

cusfbamboo.isen.p0(p, M, gamma)

Get stagnation pressure from static pressure and Mach number [1]

Parameters
  • p (float) – Pressure (Pa)

  • M (float) – Mach number

  • gamma (float) – Ratio of specific heats cp/cv

Returns

Stagnation pressure (Pa)

Return type

float

cusfbamboo.isen.T0(T, M, gamma)

Get the stangation temperature from the static temperature and Mach number [1]

Parameters
  • T (float) – Temperature (K)

  • M (float) – Mach number

  • gamma (float) – Ratio of specific heats cp/cv

Returns

Stagnation temperature (K)

Return type

float

cusfbamboo.isen.Tr(T, M, gamma, r)

Get the recovery temperature, which is the temperature you should use for convective heat transfer (i.e. q = h(Tr - Tw)). [1]

Parameters
  • T (float) – Freestream static tempreature (K)

  • M (float) – Mach number

  • gamma (float) – Ratio of specific heats cp/cv

  • r (float) – ‘Recovery factor’, usually equal to Pr^(1/3) for ‘turbulent free boundary layers’ [1]

Returns

Recovery temperature

Return type

float

cusfbamboo.isen.M_from_p(p, p0, gamma)

Mach number from static pressure and stagnation pressure.

Parameters
  • p (float) – Static pressure (Pa)

  • p0 (float) – Stagnation pressure (Pa)

  • gamma (float) – Ratio of specific heats cp/cv

Returns

Mach number

Return type

float

cusfbamboo.isen.M_from_A_subsonic(A, A_t, gamma)

Get the Mach number from the local flow area, assuming subsonic flow.

Parameters
  • A (float) – Local area (m2)

  • A_t (float) – Throat area (m2)

  • gamma (float) – Ratio of specific heats cp/cv

cusfbamboo.isen.M_from_A_supersonic(A, A_t, gamma)

Get the Mach number from the local flow area, assuming supersonic flow.

Parameters
  • A (float) – Local area (m2)

  • A_t (float) – Throat area (m2)

  • gamma (float) – Ratio of specific heats cp/cv

cusfbamboo.isen.T(T0, M, gamma)

Get local temperature from the Mach number and stagnation temperature.

Parameters
  • T0 (float) – Stagnation temperature (K)

  • M (float) – Local Mach number

  • gamma (float) – Ratio of specific heats cp/cv

Returns

Temperature (K)

Return type

float

cusfbamboo.isen.p(p0, M, gamma)

Get local pressure from the Mach number and stagnation pressure.

Parameters
  • p0 (float) – Stagnation pressure (Pa)

  • M (float) – Local Mach number

  • gamma (float) – Ratio of specific heats cp/cv

Returns

Pressure (Pa)

Return type

float

cusfbamboo.isen.get_choked_mdot(perfect_gas, chamber_conditions, At)

Get the mass flow rate through a choked nozzle.

Parameters
  • perfect_gas (PerfectGas) – Exhaust gas leaving the combustion chamber.

  • chamber_conditions (ChamberConditions) – ChamberConditions object to represent combustion chamber conditions.

  • At (float) – Throat area (m2)

cusfbamboo.isen.get_throat_area(perfect_gas, chamber_conditions, mdot)

Get the nozzle throat area, given the gas properties and combustion chamber conditions. Assumes perfect gas with isentropic flow.

Parameters
  • perfect_gas (PerfectGas) – Exhaust gas leaving the combustion chamber.

  • chamber_conditions (ChamberConditions) – ChamberConditions object to represent combustion chamber conditions.

  • mdot (float) – Mass flow rate of exhaust gas (kg/s)

Returns

Throat area (m^2)

Return type

float

cusfbamboo.isen.get_exit_area(perfect_gas, chamber_conditions, p_e, mdot)

Get the nozzle exit area, given the gas properties and combustion chamber conditions. Assumes perfect gas with isentropic flow.

Parameters
  • perfect_gas (PerfectGas) – Gas object.

  • chamber_conditions (ChamberConditions) – ChamberConditions object to represent combustion chamber conditions.

  • p_e (float) – Exit pressure (Pa)

Returns

Optimum nozzle exit area (Pa)

Return type

float

cusfbamboo.hx module

General solver for coflow and counterflow heat exchangers, using a 1-D thermal resistance model.

Notation:
  • ‘c’: Cold side (usually coolant)

  • ‘h’: Hot side (usually exhaust gas)

  • ‘w’: At the wall (e.g. T_cw is the wall temperature on the cold side)

class cusfbamboo.hx.HXSolver(T_c_in, T_h, p_c_in, cp_c, mdot_c, V_c, A_c, Rdx, extra_dQ_dx, dp_dx_f, x_start, dx, x_end)

Bases: object

Class for solving heat exchanger problems.

Parameters
  • T_c_in (float) – Coolant inlet static temperature (K)

  • T_h (callable) – Exhaust gas recovery temperature (K). Must be a function of ‘state’.

  • p_c_in (float) – Coolant inlet static pressure (Pa)

  • cp_c (callable) – Coolant isobaric specific heat capacity (J/kg/K). Must be a function of ‘state’.

  • mdot_c (float) – Coolant mass flow rate (kg/s)

  • V_c (callable) – Coolant velocity (m/s). Must be a function of ‘state’.

  • A_c (callable) – Coolant flow area (m2). Must be a function of ‘state’.

  • Rdx (callable) – List of thermal resistances [R1, R2 … etc], in the order T_cold –> T_hot. Note they need to be 1D resistances, so Qdot is per unit length. Must be a function of ‘state’.

  • extra_dQ_dx (callable) – Extra heat transfer rate (positive into the coolant), to add on (W), to represent things like fins protruding into the coolant flow. Must be a function of ‘state’.

  • dp_dx_f (callable) – Frictional pressure drop per unit length (Pa/m)

  • x_start (float) – Initial value of x to start at (m)

  • dx (float) – dx to move by for each step, corresponding to the direction that coolant flows in. Usually negative for counterflow heat exchanger (m)

  • x_end (float) – Value of x to stop at (m)

reset()

Reset our ‘state’ list to the initial conditions and set self.i to zero.

iterate()

Iterate one step at the current ‘x’ position.

step()

Move ‘dx’ onto the next x position. Make an initial guess for T_c[i+2] and p_c[i+2] based on T_c[i+1] and p_c[i+1].

run(iter_start=5, iter_each=2)

Run the simulation until we reach x >= x_end.

Parameters
  • iter_start (int, optional) – Number of iterations to use on the first gridpoint. Defaults to 5.

  • iter_each (int, optional) – Number of iterations to use on each intermediate grid point. Defaults to 1.

cusfbamboo.circuit module

Classes and functions related to thermal circuit calculations.

References

cusfbamboo.circuit.h_gas_bartz(D, cp_inf, mu_inf, Pr_inf, rho_inf, v_inf, rho_am, mu_am, mu0)

Bartz equation, using Equation (8-23) from page 312 of RPE 7th edition (Reference [1]). ‘am’ refers to the gas being at the ‘arithmetic mean’ of the wall and freestream temperatures.

Parameters
  • D (float) – Gas flow diameter (m)

  • cp_inf (float) – Specific heat capacity at constant pressure for the gas, in the freestream

  • mu_inf (float) – Absolute viscosity in the freestream

  • Pr_inf (float) – Prandtl number in the freestream

  • rho_inf (float) – Density of the gas in the freestream

  • v_inf (float) – Velocity of the gas in in the freestream

  • rho_am (float) – Density of the gas, at T = (T_wall + T_freestream)/2

  • mu_am (float) – Absolute viscosity of the gas, at T = (T_wall + T_freestream)/2

  • mu0 (float) – Absolute viscosity of the gas under stagnation conditions.

Returns

Convective heat transfer coefficient (W/m2/K), h, for the exhaust gas side (where q = h(T - T_inf)).

Return type

float

cusfbamboo.circuit.h_gas_bartz_sigma(c_star, A_t, A, p_chamber, T_chamber, M, Tw, mu0, cp0, gamma, Pr0)

Bartz heat transfer equation using the sigma correlation, from Reference [3].

Parameters
  • c_star (float) – C* efficiency ( = p_c * A_t / mdot)

  • A_t (float) – Throat area (m^2)

  • A (float) – Flow area (m^2)

  • p_chamber (float) – Chamber pressure (Pa)

  • T_chamber (float) – Chamber temperature (K)

  • M (float) – Freestream Mach number

  • Tw (float) – Wall temperature (K)

  • mu0 (float) – Absolute viscosity at stagnation conditions (Pa s)

  • cp0 (float) – Gas specific heat capacity at stagnation conditions (J/kg/K)

  • gamma (float) – Gas ratio of specific heats (cp/cv)

  • Pr0 (float) – Prandtl number at stagnation conditions

Returns

Convective heat transfer coefficient (W/m2/K), h, for the exhaust gas side (where q = h(T - T_inf)).

Return type

float

cusfbamboo.circuit.h_gas_bartz_sigma_curve(c_star, A_t, A, p_chamber, T_chamber, M, Tw, mu0, cp0, gamma, Pr0, rc_t)

Bartz heat transfer equation using the sigma correlation, from Reference [3]. Takes into account the radius of curvature at the throat.

Parameters
  • c_star (float) – C* efficiency ( = p_c * A_t / mdot)

  • A_t (float) – Throat area (m^2)

  • A (float) – Flow area (m^2)

  • p_chamber (float) – Chamber pressure (Pa)

  • T_chamber (float) – Chamber temperature (K)

  • M (float) – Freestream Mach number

  • Tw (float) – Wall temperature (K)

  • mu0 (float) – Absolute viscosity at stagnation conditions (Pa s)

  • cp0 (float) – Gas specific heat capacity at stagnation conditions (J/kg/K)

  • gamma (float) – Gas ratio of specific heats (cp/cv)

  • Pr0 (float) – Prandtl number at stagnation conditions

  • rc_t (float) – Radius of curvature at the throat

Returns

Convective heat transfer coefficient (W/m2/K), h, for the exhaust gas side (where q = h(T - T_inf)).

Return type

float

cusfbamboo.circuit.h_coolant_dittus_boelter(rho, V, D, mu, Pr, k)

Dittus-Boelter equation for convective heat transfer coefficient.

Parameters
  • rho (float) – Coolant bulk density (kg/m^3).

  • V (float) – Coolant bulk velocity (m/s)

  • D (float) – Hydraulic diameter of pipe (m)

  • mu (float) – Coolant bulk viscosity (Pa s)

  • Pr (float) – Coolant bulk Prandtl number

  • k (float) – Coolant thermal conductivity

Returns

Convective heat transfer coefficient (W/m2/K)

Return type

float

cusfbamboo.circuit.h_coolant_sieder_tate(rho, V, D, mu_bulk, mu_wall, Pr, k)

Sieder-Tate equation for convective heat transfer coefficient.

Parameters
  • rho (float) – Coolant bulk density (kg/m^3)

  • V (float) – Coolant bulk velocity (m/s)

  • D (float) – Hydraulic diameter of pipe (m)

  • mu_bulk (float) – Absolute viscosity of the coolant at the bulk temperature (Pa s).

  • mu_wall (float) – Absolute viscosity of the coolant at the wall temperature (Pa s).

  • Pr (float) – Bulk Prandtl number of the coolant.

  • k (float) – Bulk thermal conductivity of the coolant.

Returns

Convective heat transfer coefficient (W/m2/K)

Return type

float

cusfbamboo.circuit.h_coolant_gnielinski(rho, V, D, mu, Pr, k, f_darcy)

Convective heat transfer coefficient for the coolant side, using Gnielinski’s correlation. Page 41 of Reference [2].

Parameters
  • rho (float) – Coolant density (kg/m3)

  • V (float) – Coolant velocity (m/s)

  • D (float) – Hydraulic diameter (m)

  • mu (float) – Coolant viscosity (Pa s)

  • Pr (float) – Prandtl number

  • k (float) – Coolant thermal conductivity

  • f_darcy (float) – Darcy friction factor for the coolant

Returns

Convective heat transfer coefficient (W/m2/K)

Return type

float

cusfbamboo.circuit.Q_fin_adiabatic(P, Ac, k, h, L, T_b, T_inf)

Get the heat transfer rate for a fin with an adiabatic tip (Reference [5])

Parameters
  • P (float) – Fin perimeter (m)

  • Ac (float) – Fin cross sectional area (m2)

  • k (float) – Fin thermal conductivity (W/m/K)

  • h (float) – Convective heat transfer coefficient at the fin surface (W/m2/K)

  • L (float) – Fin length (m)

  • T_b (float) – Fin base temperature (K)

  • T_inf (float) – Freestream temperature of the fluid the fin is submersed in (K)

Returns

Heat transfer rate out of fin (W)

Return type

float

cusfbamboo.circuit.dQ_dA_nucleate(mu_l, h_fg, rho_l, rho_v, sigma, cp_l, T_w, T_sat, C_sf, Pr_l)

Get the heat flux due to nucleate boiling. From Rohsenow’s equation [4][6].

Parameters
  • mu_l (float) – Viscosity of the liquid phase (Pa s)

  • h_fg (float) – Enthalpy between vapour and liquid phases. h_fg = h_g - h_f. (J/kg/K)

  • rho_l (float) – Density of the liquid phase (kg/m3)

  • rho_v (float) – Density of the vapour phase (kg/m3)

  • sigma (float) – Surface tension of the liquid-vapour interface (N/m)

  • cp_l (float) – Isobaric specific heat capacity of the liquid (J/kg/K)

  • T_w (float) – Wall temperature (K)

  • T_sat (float) – Saturation temperature of the fluid (K)

  • C_sf (float) – Surface-fluid coefficient. Will be different for different material + fluid combinations. Some examples are available in [4] and [6].

  • Pr_l (float) – Prandtl number of the liquid phase

Returns

Heat flux (W/m2)

Return type

float

cusfbamboo.circuit.dQ_dA_nucleate_critical(h_fg, rho_v, sigma, rho_l)

Get the critical heat flux due to nucleate boiling, i.e. the maximum heat transfer rate that is possible. From Rohsenow’s equation [4][6].

Parameters
  • h_fg (float) – Enthalpy between vapour and liquid phases. h_fg = h_g - h_f. (J/kg/K)

  • rho_v (float) – Density of the vapour phase (kg/m3)

  • sigma (float) – Surface tension of the liquid-vapour interface (N/m)

  • rho_l (float) – Density of the liquid phase (kg/m3)

Returns

Heat flux (W/m2)

Return type

float

cusfbamboo.circuit.h_coolant_stable_film(k_vf, rho_vf, rho_v, rho_l, h_fg, cp_l, dT, mu_vf, T_w, T_sat, sigma)

Convective heat transfer coefficient for the stable-film phase of boiling heat transfer [6]. The film temperature is defined as the mean of the wall and freestream temperature, i.e. 0.5 * (T_w + T_bulk)

Parameters
  • k_vf (float) – Thermal conductivity of the vapour, evaluated at the film temperature (W/m/K)

  • rho_vf (float) – Density of the vapour, evaluated at the film temperature (kg/m3)

  • rho_v (float) – Density of the vapour, evaluated at the bulk temperature? (kg/m3)

  • rho_l (float) – Density of the liquid, evaluated at the bulk temperature? (kg/m3)

  • h_fg (float) – Enthalpy between vapour and liquid phases. h_fg = h_g - h_f. (J/kg/K)

  • cp_l (float) – Isobaric specific heat capacity of the liquid (J/kg/K)

  • dT (float) – Temperature difference between the wall and bulk (T_w - T_freestream) (K)

  • mu_vf (float) – Viscosity of the vapour, evaluated at the film temperature

  • T_w (float) – Wall temperature (K)

  • T_sat (float) – Saturated vapour temperature (K)

  • sigma (float) – Surface tension of the liquid-vapour interface (N/m)

Returns

Convective heat transfer coefficient (W/m2/K)

Return type

float

class cusfbamboo.circuit.ThermalCircuit(T1, T2, R)

Bases: object

Class for solving thermal circuits. Will solve them upon initialising.

Parameters
  • T1 (float) – Temperature at start

  • T2 (float) – Temperature at end

  • R (list) – List of resistances between T1 and T2, in the order [R_touching_T1, … , R_touching_T2]

Qdot

Heat transfer rate (positive in the direction of T1 –> T2)

Type

float

T

List of temperatures in between each resistance, including T1 and T2 at either end. i.e. [T1, …, T2].

Type

numpy.ndarray

Indices and tables