Tank and Heat Exchanger Design (fluids.geometry)

This module contains functionality for calculating parameters about different geometrical forms that come up in engineering practice.

Implemented geometry objects are tanks, helical coils, cooling towers, air coolers, compact heat exchangers, and plate and frame heat exchangers.

Additional functionality for typical catalyst/adsorbent pellet shapes is also included.

For reporting bugs, adding feature requests, or submitting pull requests, please use the GitHub issue tracker or contact the author at Caleb.Andrew.Bell@gmail.com.

Main Interfaces

class fluids.geometry.TANK(D=None, L=None, horizontal=True, sideA=None, sideB=None, sideA_a=None, sideB_a=None, sideA_f=None, sideA_k=None, sideB_f=None, sideB_k=None, sideA_a_ratio=None, sideB_a_ratio=None, L_over_D=None, V=None)[source]

Bases: object

Class representing tank volumes and levels. All parameters are also attributes.

Parameters:
Dpython:float

Diameter of the cylindrical section of the tank, [m]

Lpython:float

Length of the main cylindrical section of the tank, [m]

horizontalbool, optional

Whether or not the tank is a horizontal or vertical tank

sideApython:str, optional

The left (or bottom for vertical) head of the tank’s type; one of [None, ‘conical’, ‘ellipsoidal’, ‘torispherical’, ‘guppy’, ‘spherical’, ‘same’].

sideBpython:str, optional

The right (or top for vertical) head of the tank’s type; one of [None, ‘conical’, ‘ellipsoidal’, ‘torispherical’, ‘guppy’, ‘spherical’, ‘same’].

sideA_apython:float, optional

The distance the head as specified by sideA extends down or to the left from the main cylindrical section, [m]

sideB_apython:float, optional

The distance the head as specified by sideB extends up or to the right from the main cylindrical section, [m]

sideA_fpython:float, optional

Dimensionless dish-radius parameter for side A; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

sideA_kpython:float, optional

Dimensionless knuckle-radius parameter for side A; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

sideB_fpython:float, optional

Dimensionless dish-radius parameter for side B; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

sideB_kpython:float, optional

Dimensionless knuckle-radius parameter for side B; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

sideA_a_ratiopython:float, optional

Ratio for a parameter; can be used instead of specifying an absolute value, [-]

sideB_a_ratiopython:float, optional

Ratio for a parameter; can be used instead of specifying an absolute value, [-]

L_over_Dpython:float, optional

Ratio of length over diameter, used only when D and L are both unspecified but V is, [-]

Vpython:float, optional

Volume of the tank; solved for if specified, using sideA_a_ratio/sideB_a_ratio, sideA, sideB, horizontal, and one of L_over_D, L, or D, [m^3]

Notes

For torpsherical tank heads, the following f and k parameters are used in standards. The default is ASME F&D.

f

k

2:1 semi-elliptical

0.9

0.17

ASME F&D

1

0.06

ASME 80/6

0.8

0.06

ASME 80/10 F&D

0.8

0.1

DIN 28011

1

0.1

DIN 28013

0.8

0.154

For the following cases, numerical integrals are used.

V_horiz_spherical V_horiz_torispherical SA_partial_horiz_spherical_head SA_partial_horiz_ellipsoidal_head SA_partial_horiz_guppy_head SA_partial_horiz_torispherical_head

Examples

Total volume of a tank:

>>> TANK(D=1.2, L=4, horizontal=False).V_total
4.523893421169302

Volume of a tank at a given height:

>>> TANK(D=1.2, L=4, horizontal=False).V_from_h(.5)
0.5654866776461628

Height of liquid for a given volume:

>>> TANK(D=1.2, L=4, horizontal=False).h_from_V(.5)
0.442097064

Surface area of a tank with a conical head:

>>> T1 = TANK(V=10, L_over_D=0.7, sideB='conical', sideB_a=0.5)
>>> T1.A, T1.A_sideA, T1.A_sideB, T1.A_lateral
(24.94775907, 5.118555, 5.497246, 14.331956)

Solving for tank volumes, first horizontal, then vertical:

>>> TANK(D=10., horizontal=True, sideA='conical', sideB='conical', V=500).L
4.699531
>>> TANK(L=4.69953105701, horizontal=True, sideA='conical', sideB='conical', V=500).D
9.9999999
>>> TANK(L_over_D=0.469953105701, horizontal=True, sideA='conical', sideB='conical', V=500).L
4.6995310
>>> TANK(D=10., horizontal=False, sideA='conical', sideB='conical', V=500).L
4.699531
>>> TANK(L=4.69953105701, horizontal=False, sideA='conical', sideB='conical', V=500).D
9.99999999
>>> TANK(L_over_D=0.469953105701, horizontal=False, sideA='conical', sideB='conical', V=500).L
4.699531057
Attributes:
h_maxpython:float

Height of the tank, [m]

V_totalpython:float

Total volume of the tank as calculated [m^3]

sideA_Vpython:float

Volume of only sideA [m^3]

sideB_Vpython:float

Volume of only sideB [m^3]

lateral_Vpython:float

Volume of cylindrical section of tank [m^3]

Apython:float

Total surface area of the tank, [m^2]

A_sideApython:float

Surface area of sideA, [m^2]

A_sideBpython:float

Surface area of sideB, [m^2]

A_lateralpython:float

Surface area of the lateral side, [m^2]

A_sideA_extrapython:float

Additional surface area of sideA beyond that of a flat disk, [m^2]

A_sideB_extrapython:float

Additional surface area of sideB beyond that of a flat disk, [m^2]

tablebool

Whether or not a table of heights-volumes has been generated

heightsndarray

Array of heights between 0 and h_max, [m]

volumesndarray

Array of volumes calculated from the heights, [m^3]

c_forwardndarray

Coefficients for the Chebyshev approximations in calculating V from h, [-]

c_backwardndarray

Coefficients for the Chebyshev approximations in calculating h from V, [-]

Methods

A_cross_sectional(h[, method])

Method to calculate the cross-sectional liquid surface area from which gas can evolve in a fully defined tank given a specified height h.

SA_from_h(h[, method])

Method to calculate the volume of liquid in a fully defined tank given a specified height h.

V_from_h(h[, method])

Method to calculate the volume of liquid in a fully defined tank given a specified height h.

add_thickness(thickness[, sideA_thickness, ...])

Method to create a new tank instance with the same parameters as itself, except with an added thickness to it.

from_two_specs(spec0, spec1[, spec0_name, ...])

Method to create a new tank instance according to two specifications which are not direct geometry parameters.

h_from_V(V[, method])

Method to calculate the height of liquid in a fully defined tank given a specified volume of liquid in it V.

set_chebyshev_approximators([deg_forward, ...])

Method to derive and set coefficients for chebyshev polynomial function approximation of the height-volume and volume-height relationship.

set_misc()

Set more parameters, after the tank is better defined than in the __init__ function.

set_table([n, dx])

Method to set an interpolation table of liquids levels versus volumes in the tank, for a fully defined tank.

A_cross_sectional(h, method='full')[source]

Method to calculate the cross-sectional liquid surface area from which gas can evolve in a fully defined tank given a specified height h. h must be under the maximum height. This is calculated by numeric differentiation for most cases.

Parameters:
hpython:float

Height specified, [m]

methodpython:str, optional

‘full’ (calculated rigorously) or ‘chebyshev’, [-]

Returns:
A_crosspython:float

Surface area of liquid in the tank up to the specified height, [m^2]

SA_from_h(h, method='full')[source]

Method to calculate the volume of liquid in a fully defined tank given a specified height h. h must be under the maximum height.

Parameters:
hpython:float

Height specified, [m]

methodpython:str, optional

‘full’ (calculated rigorously) ; nothing else is implemented

Returns:
SApython:float

Surface area of liquid in the tank up to the specified height, [m^2]

V_from_h(h, method='full')[source]

Method to calculate the volume of liquid in a fully defined tank given a specified height h. h must be under the maximum height. If the method is ‘chebyshev’, and the coefficients have not yet been calculated, they are created by calling set_chebyshev_approximators.

Parameters:
hpython:float

Height specified, [m]

methodpython:str

One of ‘full’ (calculated rigorously) or ‘chebyshev’

Returns:
Vpython:float

Volume of liquid in the tank up to the specified height, [m^3]

add_thickness(thickness, sideA_thickness=None, sideB_thickness=None)[source]

Method to create a new tank instance with the same parameters as itself, except with an added thickness to it. This is useful to obtain ex. the inside of a tank and the outside; their different in volumes is the volume of the shell, and could be used to determine weight.

Parameters:
thicknesspython:float

Thickness to add to the tank diameter, [m]

sideA_thicknesspython:float, optional

The thickness to add to the sideA head; if not specified, it will be thickness, [m]

sideB_thicknesspython:float, optional

The thickness to add to the sideB head; if not specified, it will be thickness, [m]

Returns:
TANKTANK

Tank object, [-]

Notes

Be careful not to specify a negative thickness larger than the heads’ lengths, or the head will become concave! The same applies to adding a thickness to convex heads - they can become convex.

chebyshev = False
static from_two_specs(spec0, spec1, spec0_name='V', spec1_name='A_cross', h=None, horizontal=True, sideA=None, sideB=None, sideA_a=None, sideB_a=None, sideA_f=None, sideA_k=None, sideB_f=None, sideB_k=None, sideA_a_ratio=None, sideB_a_ratio=None)[source]

Method to create a new tank instance according to two specifications which are not direct geometry parameters.

The allowable options are ‘V’, ‘SA’, ‘V_partial’, ‘SA_partial’, and ‘A_cross’, the later three of which require h to be specified.

Parameters:
spec0python:float

Goal for spec0_name, [-]

spec1python:float

Goal for spec1_name, [-]

spec0_namepython:str

One of ‘V’, ‘SA’, ‘V_partial’, ‘SA_partial’, and ‘A_cross’ [-]

spec1_namepython:str

One of ‘V’, ‘SA’, ‘V_partial’, ‘SA_partial’, and ‘A_cross’ [-]

hpython:float

Height at which to calculate the specs, [m]

horizontalbool, optional

Whether or not the tank is a horizontal or vertical tank

sideApython:str, optional

The left (or bottom for vertical) head of the tank’s type; one of [None, ‘conical’, ‘ellipsoidal’, ‘torispherical’, ‘guppy’, ‘spherical’, ‘same’].

sideBpython:str, optional

The right (or top for vertical) head of the tank’s type; one of [None, ‘conical’, ‘ellipsoidal’, ‘torispherical’, ‘guppy’, ‘spherical’, ‘same’].

sideA_apython:float, optional

The distance the head as specified by sideA extends down or to the left from the main cylindrical section, [m]

sideB_apython:float, optional

The distance the head as specified by sideB extends up or to the right from the main cylindrical section, [m]

sideA_fpython:float, optional

Dimensionless dish-radius parameter for side A; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

sideA_kpython:float, optional

Dimensionless knuckle-radius parameter for side A; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

sideB_fpython:float, optional

Dimensionless dish-radius parameter for side B; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

sideB_kpython:float, optional

Dimensionless knuckle-radius parameter for side B; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

Returns:
TANKTANK

Tank object at solved specifications, [-]

Notes

Limited testing has been done on this method. The bounds are D between 0.1 mm and 10 km, with L_D ratios of 1e-4 to 1e4.

h_from_V(V, method='spline')[source]

Method to calculate the height of liquid in a fully defined tank given a specified volume of liquid in it V. V must be under the maximum volume. If the method is ‘spline’, and the interpolation table is not yet defined, creates it by calling the method set_table. If the method is ‘chebyshev’, and the coefficients have not yet been calculated, they are created by calling set_chebyshev_approximators.

Parameters:
Vpython:float

Volume of liquid in the tank up to the desired height, [m^3]

methodpython:str

One of ‘spline’, ‘chebyshev’, or ‘brenth’

Returns:
hpython:float

Height of liquid at which the volume is as desired, [m]

model_inputs = ('D', 'L', 'horizontal', 'sideA', 'sideB', 'sideA_a', 'sideB_a', 'sideA_f', 'sideA_k', 'sideB_f', 'sideB_k', 'sideA_a_ratio', 'sideB_a_ratio', 'L_over_D', 'V')
set_chebyshev_approximators(deg_forward=50, deg_backwards=200)[source]

Method to derive and set coefficients for chebyshev polynomial function approximation of the height-volume and volume-height relationship.

A single set of chebyshev coefficients is used for the entire height- volume and volume-height relationships respectively.

The forward relationship, V_from_h, requires far fewer coefficients in its fit than the reverse to obtain the same relative accuracy.

Optionally, deg_forward or deg_backwards can be set to None to try to automatically fit the series to machine precision.

Parameters:
deg_forwardpython:int, optional

The degree of the chebyshev polynomial to be created for the V_from_h curve, [-]

deg_backwardspython:int, optional

The degree of the chebyshev polynomial to be created for the h_from_V curve, [-]

set_misc()[source]

Set more parameters, after the tank is better defined than in the __init__ function.

Notes

Two of D, L, and L_over_D must be known when this function runs. The other one is set from the other two first thing in this function. a_ratio parameters are used to calculate a values for the heads here, if applicable. Radius is calculated here. Maximum tank height is calculated here. V_total is calculated here.

set_table(n=100, dx=None)[source]

Method to set an interpolation table of liquids levels versus volumes in the tank, for a fully defined tank. Normally run by the h_from_V method, this may be run prior to its use with a custom specification. Either the number of points on the table, or the vertical distance between steps may be specified.

Parameters:
npython:float, optional

Number of points in the interpolation table, [-]

dxpython:float, optional

Vertical distance between steps in the interpolation table, [m]

table = False
fluids.geometry.V_tank(D, L, horizontal=True, sideA=None, sideB=None, sideA_a=0.0, sideB_a=0.0, sideA_f=None, sideA_k=None, sideB_f=None, sideB_k=None)[source]

Calculates the total volume of a vertical or horizontal tank with different head types.

Parameters:
Dpython:float

Diameter of the cylindrical section of the tank, [m]

Lpython:float

Length of the main cylindrical section of the tank, [m]

horizontalbool, optional

Whether or not the tank is a horizontal or vertical tank

sideApython:str, optional

The left (or bottom for vertical) head of the tank’s type; one of [None, ‘conical’, ‘ellipsoidal’, ‘torispherical’, ‘guppy’, ‘spherical’].

sideBpython:str, optional

The right (or top for vertical) head of the tank’s type; one of [None, ‘conical’, ‘ellipsoidal’, ‘torispherical’, ‘guppy’, ‘spherical’].

sideA_apython:float, optional

The distance the head as specified by sideA extends down or to the left from the main cylindrical section, [m]

sideB_apython:float, optional

The distance the head as specified by sideB extends up or to the right from the main cylindrical section, [m]

sideA_fpython:float, optional

Dimensionless dish-radius parameter for side A; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

sideA_kpython:float, optional

Dimensionless knuckle-radius parameter for side A; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

sideB_fpython:float, optional

Dimensionless dish-radius parameter for side B; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

sideB_kpython:float, optional

Dimensionless knuckle-radius parameter for side B; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

Returns:
Vpython:float

Total volume [m^3]

sideA_Vpython:float

Volume of only sideA [m^3]

sideB_Vpython:float

Volume of only sideB [m^3]

lateral_Vpython:float

Volume of cylindrical section of tank [m^3]

Examples

>>> V_tank(D=1.5, L=5., horizontal=False, sideA='conical',
... sideB='conical', sideA_a=2., sideB_a=1.)
(10.602875205865551, 1.1780972450961726, 0.5890486225480863, 8.835729338221293)
fluids.geometry.V_from_h(h, D, L, horizontal=True, sideA=None, sideB=None, sideA_a=0, sideB_a=0, sideA_f=None, sideA_k=None, sideB_f=None, sideB_k=None)[source]

Calculates partially full volume of a vertical or horizontal tank with different head types according to [1].

If the height specified is above the height of the tank, it is truncated to the top of the tank.

Parameters:
hpython:float

Height of the liquid in the tank, [m]

Dpython:float

Diameter of the cylindrical section of the tank, [m]

Lpython:float

Length of the main cylindrical section of the tank, [m]

horizontalbool, optional

Whether or not the tank is a horizontal or vertical tank

sideApython:str, optional

The left (or bottom for vertical) head of the tank’s type; one of [None, ‘conical’, ‘ellipsoidal’, ‘torispherical’, ‘guppy’, ‘spherical’].

sideBpython:str, optional

The right (or top for vertical) head of the tank’s type; one of [None, ‘conical’, ‘ellipsoidal’, ‘torispherical’, ‘guppy’, ‘spherical’].

sideA_apython:float, optional

The distance the head as specified by sideA extends down or to the left from the main cylindrical section, [m]

sideB_apython:float, optional

The distance the head as specified by sideB extends up or to the right from the main cylindrical section, [m]

sideA_fpython:float, optional

Dimensionless dish-radius parameter for side A; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

sideA_kpython:float, optional

Dimensionless knuckle-radius parameter for side A; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

sideB_fpython:float, optional

Dimensionless dish-radius parameter for side B; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

sideB_kpython:float, optional

Dimensionless knuckle-radius parameter for side B; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

Returns:
Vpython:float

Volume up to h [m^3]

References

[1]

Jones, D. “Compute Fluid Volumes in Vertical Tanks.” Chemical Processing. December 18, 2003. http://www.chemicalprocessing.com/articles/2003/193/

Examples

>>> V_from_h(h=7, D=1.5, L=5., horizontal=False, sideA='conical',
... sideB='conical', sideA_a=2., sideB_a=1.)
10.013826583317465
fluids.geometry.SA_tank(D, L, sideA=None, sideB=None, sideA_a=0, sideB_a=0, sideA_f=None, sideA_k=None, sideB_f=None, sideB_k=None)[source]

Calculates the surface are of a cylindrical tank with optional heads. In the degenerate case of being provided with only D and L, provides the surface area of a cylinder.

Parameters:
Dpython:float

Diameter of the cylindrical section of the tank, [m]

Lpython:float

Length of the main cylindrical section of the tank, [m]

sideApython:str, optional

The left (or bottom for vertical) head of the tank’s type; one of [None, ‘conical’, ‘ellipsoidal’, ‘torispherical’, ‘guppy’, ‘spherical’].

sideBpython:str, optional

The right (or top for vertical) head of the tank’s type; one of [None, ‘conical’, ‘ellipsoidal’, ‘torispherical’, ‘guppy’, ‘spherical’].

sideA_apython:float, optional

The distance the head as specified by sideA extends down or to the left from the main cylindrical section, [m]

sideB_apython:float, optional

The distance the head as specified by sideB extends up or to the right from the main cylindrical section, [m]

sideA_fpython:float, optional

Dish-radius parameter for side A; fD = dish radius [1/m]

sideA_kpython:float, optional

knuckle-radius parameter for side A; kD = knuckle radius [1/m]

sideB_fpython:float, optional

Dish-radius parameter for side B; fD = dish radius [1/m]

sideB_kpython:float, optional

knuckle-radius parameter for side B; kD = knuckle radius [1/m]

Returns:
SApython:float

Surface area of the tank [m^2]

sideA_SApython:float

Surface area of only sideA [m^2]

sideB_SApython:float

Surface area of only sideB [m^2]

lateral_SApython:float

Surface area of cylindrical section of tank [m^2]

Examples

Cylinder, Spheroid, Long Cones, and spheres. All checked.

>>> SA_tank(D=2, L=2)[0]
18.84955592153876
>>> SA_tank(D=1., L=0, sideA='ellipsoidal', sideA_a=2, sideB='ellipsoidal',
... sideB_a=2)[0]
10.124375616183
>>> SA_tank(D=1., L=5, sideA='conical', sideA_a=2, sideB='conical',
... sideB_a=2)[0]
22.18452243965
>>> SA_tank(D=1., L=5, sideA='spherical', sideA_a=0.5, sideB='spherical',
... sideB_a=0.5)[0]
18.8495559215
fluids.geometry.SA_from_h(h, D, L, horizontal=True, sideA=None, sideB=None, sideA_a=0.0, sideB_a=0.0, sideA_f=None, sideA_k=None, sideB_f=None, sideB_k=None)[source]

Calculates partially full wetted surface area of a vertical or horizontal tank with different head types according to [1].

Parameters:
hpython:float

Height of the liquid in the tank, [m]

Dpython:float

Diameter of the cylindrical section of the tank, [m]

Lpython:float

Length of the main cylindrical section of the tank, [m]

horizontalbool, optional

Whether or not the tank is a horizontal or vertical tank

sideApython:str, optional

The left (or bottom for vertical) head of the tank’s type; one of [None, ‘conical’, ‘ellipsoidal’, ‘torispherical’, ‘guppy’, ‘spherical’].

sideBpython:str, optional

The right (or top for vertical) head of the tank’s type; one of [None, ‘conical’, ‘ellipsoidal’, ‘torispherical’, ‘guppy’, ‘spherical’].

sideA_apython:float, optional

The distance the head as specified by sideA extends down or to the left from the main cylindrical section, [m]

sideB_apython:float, optional

The distance the head as specified by sideB extends up or to the right from the main cylindrical section, [m]

sideA_fpython:float, optional

Dimensionless dish-radius parameter for side A; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

sideA_kpython:float, optional

Dimensionless knuckle-radius parameter for side A; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

sideB_fpython:float, optional

Dimensionless dish-radius parameter for side B; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

sideB_kpython:float, optional

Dimensionless knuckle-radius parameter for side B; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

Returns:
SApython:float

Wetted wall surface area up to h [m^3]

References

[1]

Jones, D. “Calculating Tank Wetted Area.” Text. Chemical Processing. April 2017. https://www.chemicalprocessing.com/assets/Uploads/calculating-tank-wetted-area.pdf http://www.chemicalprocessing.com/articles/2003/193/

Examples

>>> SA_from_h(h=7, D=1.5, L=5., horizontal=False, sideA='conical',
... sideB='conical', sideA_a=2., sideB_a=1.)
28.59477853914843
class fluids.geometry.HelicalCoil(Dt, Do=None, pitch=None, H=None, N=None, H_total=None, Do_total=None, Di=None)[source]

Bases: object

Class representing a helical coiled tube, as are found in many heated tanks and some small nuclear reactors. All parameters are also attributes.

One set of the following parameters is required; inner tube diameter is optional.

  • Tube outer diameter, coil outer diameter, pitch, number of coil turns

  • Tube outer diameter, coil outer diameter, pitch, height

  • Tube outer diameter, coil outer diameter, number of coil turns, height

Parameters:
Dtpython:float

Outer diameter of the tube wound to make up the helical spiral, [m]

Dopython:float

Diameter of the spiral as measured from the center of the coil on one side to the center of the coil on the other side, [m]

Do_totalpython:float, optional

Diameter of the spiral as measured from one edge of the tube to the other edge; equal to Do + Dt; either Do or Do_total may be specified and the other will be calculated [m]

pitchpython:float, optional

Height change from one coil to the next as measured from the middles of the tube, [m]

Hpython:float, optional

Height of the spiral, as measured from the middle of the bottom of the tube to the middle of the top of the tube, [m]

H_totalpython:float, optional

Height of the spiral as measured from one edge of the tube to the other edge; equal to H_total + Dt; either may be specified and the other will be calculated [m]

Npython:float, optional

Number of coil turns; may be specified along with pitch instead of specifying H or H_total, [-]

Dipython:float, optional

Inner diameter of the tube; if specified, inside and annulus properties will be calculated, [m]

Notes

Do must be larger than Dt.

References

[1]

El-Genk, Mohamed S., and Timothy M. Schriener. “A Review and Correlations for Convection Heat Transfer and Pressure Losses in Toroidal and Helically Coiled Tubes.” Heat Transfer Engineering 0, no. 0 (June 7, 2016): 1-28. doi:10.1080/01457632.2016.1194693.

Examples

>>> C1 = HelicalCoil(Do=30, H=20, pitch=5, Dt=2)
>>> C1.N, C1.tube_length, C1.surface_area
(4.0, 377.5212621504738, 2372.0360474917497)

Same coil, with the inputs one would physically measure from the coil, and a specified inlet diameter:

>>> C1 = HelicalCoil(Do_total=32, H_total=22, pitch=5, Dt=2, Di=1.8)
>>> C1.N, C1.tube_length, C1.surface_area
(4.0, 377.5212621504738, 2372.0360474917497)
>>> C1.inner_surface_area, C1.inlet_area, C1.inner_volume, C1.total_volume, C1.annulus_volume
(2134.832442742575, 2.5446900494077327, 960.6745992341587, 1186.0180237458749, 225.3434245117162)
Attributes:
tube_circumferencepython:float

Circumference of the tube as measured though its center, not inner or outer edges; \(C = \pi D_o\), [m]

tube_lengthpython:float

Length of tube used to make the helical coil; \(L = \sqrt{(\pi D_o\cdot N)^2 + H^2}\), [m]

surface_areapython:float

Surface area of the outer surface of the helical coil; \(A_t = \pi D_t L\), [m^2]

inner_surface_areapython:float

Surface area of the inner surface of the helical coil; calculated if Di is supplied; \(A_{inside} = \pi D_i L\), [m^2]

inlet_areapython:float

Area of the inlet to the helical coil; calculated if Di is supplied; \(A_{inlet} = \frac{\pi}{4} D_i^2\), [m^2]

inner_volumepython:float

Volume of the tube as would be filled by a fluid, useful for weight calculations; calculated if Di is supplied; \(V_{inside} = A_i L\), [m^3]

annulus_areapython:float

Area of the annulus (wall of the pipe); calculated if Di is supplied; \(A_a = \frac{\pi}{4} (D_t^2 - D_i^2)\), [m^2]

annulus_volumepython:float

Volume of the annulus (wall of the pipe); calculated if Di is supplied, useful for weight calculations; \(V_a = A_a L\), [m^3]

total_volumepython:float

Total volume occupied by the pipe and the fluid inside it; \(V = D_t L\), [m^3]

helix_anglepython:float

Angle between the pitch and coil diameter; used in some calculations; \(\alpha = \arctan \left(\frac{p_t}{\pi D_o}\right)\), [radians]

curvaturepython:float

Coil curvature, useful in some calculations; \(\delta = \frac{D_t}{D_o[1 + 4\pi^2 \tan^2(\alpha)]}\), [-]

class fluids.geometry.PlateExchanger(amplitude, wavelength, chevron_angle=45, chevron_angles=None, width=None, length=None, thickness=None, d_port=None, plates=None)[source]

Bases: object

Class representing a plate heat exchanger with sinusoidal ridges. All parameters are also attributes.

Parameters:
amplitudepython:float

Half the height of the wave of the ridges, [m]

wavelengthpython:float

Distance between the bottoms of two of the ridges (sometimes called pitch), [m]

chevron_anglepython:float, optional

Angle of the plate corrugations with respect to the vertical axis (the direction of flow if the plates were straight), between 0 and 90, [degrees]

chevron_anglespython:tuple(2), optional

Many plate exchangers use two alternating patterns; for those cases provide tuple of the two angles for that situation and the argument chevron_angle is ignored, [degrees]

widthpython:float, optional

Width of the plates in the heat exchanger, between the gaskets, [m]

lengthpython:float, optional

Length of the heat exchanger as measured from one port to the other, excluding the diameter of the ports themselves (little useful heat transfer happens there), [m]

thicknesspython:float, optional

Thickness of the metal making up the plates, [m]

d_portpython:float, optional

The diameter of the ports in the plates, [m]

platespython:int, optional

The number of plates in the heat exchanger, including the two not used for heat transfer at the beginning and end [-]

Notes

Only wavelength and amplitude are required as inputs to this function.

References

[1]

Amalfi, Raffaele L., Farzad Vakili-Farahani, and John R. Thome. “Flow Boiling and Frictional Pressure Gradients in Plate Heat Exchangers. Part 1: Review and Experimental Database.” International Journal of Refrigeration 61 (January 2016): 166-84. doi:10.1016/j.ijrefrig.2015.07.010.

Examples

>>> PlateExchanger(amplitude=5E-4, wavelength=3.7E-3, length=1.2, width=.3,
... d_port=.05, plates=51)
<Plate heat exchanger, amplitude=0.0005 m, wavelength=0.0037 m, chevron_angles=45/45 degrees, area enhancement factor=1.16119, width=0.3 m, length=1.2 m, port diameter=0.05 m, heat transfer area=20.4833 m^2, 51 plates>
Attributes:
chevron_anglespython:tuple(2)

The two specified angles (repeated value if only one specified), [degrees]

chevron_anglepython:float

The averaged angle of the chevrons, [degrees]

inclination_anglepython:float

90 - chevron_angle, used in many publications instead of chevron_angle, [degrees]

plate_corrugation_aspect_ratiopython:float

The aspect ratio of the corrugations \(\gamma = \frac{4a}{\lambda}\), [-]

plate_enlargement_factorpython:float

The extra surface area multiplier as compared to a flat plate caused the corrugations, [-]

D_eqpython:float

Equivalent diameter of the channels, \(D_{eq} = 4a\) [m]

D_hydraulicpython:float

Hydraulic diameter of the channels, \(D_{hyd} = \frac{4a}{\phi}\) [m]

length_portpython:float

Port center to port center along the direction of flow, [m]

A_plate_surfacepython:float

The surface area of one plate in the heat exchanger, including the extra due to corrugations (excluding the bit between the ports), \(A_p = L\cdot W\cdot \phi\) [m^2]

A_heat_transferpython:float

The total surface area available for heat transfer in the exchanger, the multiple of A_plate_surface by the number of plates after removing the two on the edges, [m^2]

A_channel_flowpython:float

The area for the fluid to flow in one channel, \(W\cdot b\) [m^2]

channelspython:int

The number of plates minus one, [-]

channels_per_fluidpython:int

Half the number of total channels, [-]

property plate_exchanger_identifier

Method to create an identifying string in format ‘L’ + wavelength + ‘A’ + amplitude + ‘B’ + chevron angle-chevron angle.

Wavelength and amplitude are specified in units of mm and rounded to two decimal places.

class fluids.geometry.AirCooledExchanger(tube_rows, tube_passes, tubes_per_row, tube_length, tube_diameter, fin_thickness, angle=None, pitch=None, pitch_parallel=None, pitch_normal=None, pitch_ratio=None, fin_diameter=None, fin_height=None, fin_density=None, fin_interval=None, parallel_bays=1, bundles_per_bay=1, fans_per_bay=1, corbels=False, tube_thickness=None, fan_diameter=None)[source]

Bases: object

Class representing the geometry of an air cooled heat exchanger with one or more tube bays, fans, or bundles. All parameters are also attributes.

The minimum information required to describe an air cooler is as follows:

  • tube_rows

  • tube_passes

  • tubes_per_row

  • tube_length

  • tube_diameter

  • fin_thickness

Two of angle, pitch, pitch_parallel, and pitch_normal (pitch_ratio may take the place of pitch)

Either fin_diameter or fin_height. Either fin_density or fin_interval.

Parameters:
tube_rowspython:int

Number of tube rows per bundle, [-]

tube_passespython:int

Number of tube passes (times the fluid travels across one tube length), [-]

tubes_per_rowpython:float

Number of tubes per row per bundle, [-]

tube_lengthpython:float

Total length of the tube bundle tubes, [m]

tube_diameterpython:float

Diameter of the bare tube, [m]

fin_thicknesspython:float

Thickness of the fins, [m]

anglepython:float, optional

Angle of the tube layout, [degrees]

pitchpython:float, optional

Shortest distance between tube centers; defined in relation to the flow direction only, [m]

pitch_parallelpython:float, optional

Distance between tube center along a line parallel to the flow; has been called longitudinal pitch, pp, s2, SL, and p2, [m]

pitch_normalpython:float, optional

Distance between tube centers in a line 90° to the line of flow; has been called the transverse pitch, pn, s1, ST, and p1, [m]

pitch_ratiopython:float, optional

Ratio of the pitch to bare tube diameter, [-]

fin_diameterpython:float, optional

Outer diameter of each tube after including the fin on both sides, [m]

fin_heightpython:float, optional

Height above bare tube of the tube fins, [m]

fin_densitypython:float, optional

Number of fins per meter of tube, [1/m]

fin_intervalpython:float, optional

Space between each fin, including the thickness of one fin at its base, [m]

parallel_bayspython:int, optional

Number of bays in the unit, [-]

bundles_per_baypython:int, optional

Number of tube bundles per bay, [-]

fans_per_baypython:int, optional

Number of fans per bay, [-]

corbelsbool, optional

Whether or not the air cooler has corbels, which increase the air velocity by adding half a tube to the sides for the case of non-rectangular tube layouts, [-]

tube_thicknesspython:float, optional

Thickness of the bare metal tubes, [m]

fan_diameterpython:float, optional

Diameter of air cooler fan, [m]

References

[1]

Schlunder, Ernst U, and International Center for Heat and Mass Transfer. Heat Exchanger Design Handbook. Washington: Hemisphere Pub. Corp., 1983.

Examples

>>> from scipy.constants import inch
>>> AC = AirCooledExchanger(tube_rows=4, tube_passes=4, tubes_per_row=56, tube_length=10.9728,
... tube_diameter=1*inch, fin_thickness=0.013*inch, fin_density=10/inch,
... angle=30, pitch=2.5*inch, fin_height=0.625*inch, tube_thickness=0.00338,
... bundles_per_bay=2, parallel_bays=3, corbels=True)
Attributes:
bare_lengthpython:float

Length of bare tube between two fins \(\text{bare length} = \text{fin interval} - t_{fin}\), [m]

tubes_per_bundlepython:float

Total number of tubes per bundle \(N_{tubes/bundle} = N_{tubes/row} \cdot N_{rows}\), [-]

tubes_per_baypython:float

Total number of tubes per bay \(N_{tubes/bay} = N_{tubes/bundle} \cdot N_{bundles/bay}\), [-]

tubespython:float

Total number of tubes in all bundles in all bays combined \(N_{tubes} = N_{tubes/bay} \cdot N_{bays}\), [-]

pitch_diagonalpython:float

Distance between tube centers in a diagonal line between one normal tube and one parallel tube; \(s_D = \left[s_L^2 + \left(\frac{s_T}{2}\right)^2\right]^{0.5}\), [m]

A_bare_tube_per_tubepython:float

Area of the bare tube including the portion hidden by the fin per tube \(A_{bare,total/tube} = \pi D_{tube} L_{tube}\), [m^2]

A_bare_tube_per_rowpython:float

Area of the bare tube including the portion hidden by the fin per tube row \(A_{bare,total/row} = \pi D_{tube} L_{tube} N_{tubes/row}\), [m^2]

A_bare_tube_per_bundlepython:float

Area of the bare tube including the portion hidden by the fin per bundle \(A_{bare,total/bundle} = \pi D_{tube} L_{tube} N_{tubes/bundle}\), [m^2]

A_bare_tube_per_baypython:float

Area of the bare tube including the portion hidden by the fin per bay \(A_{bare,total/bay} = \pi D_{tube} L_{tube} N_{tubes/bay}\), [m^2]

A_bare_tubepython:float

Area of the bare tube including the portion hidden by the fin per in all bundles and bays combined \(A_{bare,total} = \pi D_{tube} L_{tube} N_{tubes}\), [m^2]

A_tube_showing_per_tubepython:float

Area of the bare tube which is exposed per tube \(A_{bare, showing/tube} = \pi D_{tube} L_{tube} \left(1 - \frac{t_{fin}} {\text{fin interval}} \right)\), [m^2]

A_tube_showing_per_rowpython:float

Area of the bare tube which is exposed per tube row, [m^2]

A_tube_showing_per_bundlepython:float

Area of the bare tube which is exposed per bundle, [m^2]

A_tube_showing_per_baypython:float

Area of the bare tube which is exposed per bay, [m^2]

A_tube_showingpython:float

Area of the bare tube which is exposed in all bundles and bays combined, [m^2]

A_per_finpython:float

Surface area per fin \(A_{fin} = 2 \frac{\pi}{4} (D_{fin}^2 - D_{tube}^2) + \pi D_{fin} t_{fin}\), [m^2]

A_fin_per_tubepython:float

Surface area of all fins per tube \(A_{fin/tube} = N_{fins/m} L_{tube} A_{fin}\), [m^2]

A_fin_per_rowpython:float

Surface area of all fins per row, [m^2]

A_fin_per_bundlepython:float

Surface area of all fins per bundle, [m^2]

A_fin_per_baypython:float

Surface area of all fins per bay, [m^2]

A_finpython:float

Surface area of all fins in all bundles and bays combined, [m^2]

A_per_tubepython:float

Surface area of combined finned and non-fined area exposed for heat transfer per tube \(A_{tube} = A_{bare, showing/tube} + A_{fin/tube}\), [m^2]

A_per_rowpython:float

Surface area of combined finned and non-finned area exposed for heat transfer per tube row, [m^2]

A_per_bundlepython:float

Surface area of combined finned and non-finned area exposed for heat transfer per tube bundle, [m^2]

A_per_baypython:float

Surface area of combined finned and non-finned area exposed for heat transfer per bay, [m^2]

Apython:float

Surface area of combined finned and non-finned area exposed for heat transfer in all bundles and bays combined, [m^2]

A_increasepython:float

Ratio of actual surface area to bare tube surface area \(A_{increase} = \frac{A_{tube}}{A_{bare, total/tube}}\), [-]

A_tube_flowpython:float

The area for the fluid to flow in one tube, \(\pi/4\cdot D_i^2\), [m^2]

A_tube_flow_per_rowpython:float

The area for the fluid to flow in one row, \(\pi/4\cdot D_i^2 N_{tubes/row}\), [m^2]

A_tube_flow_per_bundlepython:float

The area for the fluid to flow in one bundle, \(\pi/4\cdot D_i^2 N_{tubes/bundle}\), [m^2]

A_tube_flow_per_baypython:float

The area for the fluid to flow in one bay, \(\pi/4\cdot D_i^2 N_{tubes/bay}\), [m^2]

A_tube_flow_totalpython:float

The area for the fluid to flow in all tubes combined, as if there were a single pass [m^2]

channelspython:int

The number of tubes the fluid flows through at the inlet header, [-]

tube_volume_per_tubepython:float

Fluid volume per tube inside \(V_{tube, flow} = \frac{\pi}{4} D_{i}^2 L_{tube}\), [m^3]

tube_volume_per_rowpython:float

Fluid volume of tubes per row, [m^3]

tube_volume_per_bundlepython:float

Fluid volume of tubes per bundle, [m^3]

tube_volume_per_baypython:float

Fluid volume of tubes per bay, [m^3]

tube_volumepython:float

Fluid volume of tubes in all bundles and bays combined, [m^3]

A_diagonal_per_bundlepython:float

Air flow area along the diagonal plane per bundle \(A_d = 2 N_{tubes/row} L_{tube} (P_d - D_{tube} - 2 N_{fins/m} h_{fin} t_{fin}) + A_\text{extra,side}\), [m^2]

A_normal_per_bundlepython:float

Air flow area along the normal (transverse) plane; this is normally the minimum flow area, except for some staggered configurations \(A_t = N_{tubes/row} L_{tube} (P_t - D_{tube} - 2 N_{fins/m} h_{fin} t_{fin}) + A_\text{extra,side}\), [m^2]

A_min_per_bundlepython:float

Minimum air flow area per bundle; this is the characteristic area for velocity calculation in most finned tube convection correlations \(A_{min} = min(A_d, A_t)\), [m^2]

A_min_per_baypython:float

Minimum air flow area per bay, [m^2]

A_minpython:float

Minimum air flow area, [m^2]

A_face_per_bundlepython:float

Face area per bundle \(A_{face} = P_{T} (1+N_{tubes/row}) L_{tube}\); if corbels are used, add 0.5 to tubes/row instead of 1, [m^2]

A_face_per_baypython:float

Face area per bay, [m^2]

A_facepython:float

Total face area, [m^2]

flow_area_contraction_ratiopython:float

Ratio of A_min to A_face, [-]

model_inputs = ('tube_rows', 'tube_passes', 'tubes_per_row', 'tube_length', 'tube_diameter', 'fin_thickness', 'angle', 'pitch', 'pitch_parallel', 'pitch_normal', 'fin_diameter', 'fin_height', 'fin_density', 'fin_interval', 'parallel_bays', 'bundles_per_bay', 'fans_per_bay', 'corbels', 'tube_thickness', 'fan_diameter')
class fluids.geometry.HyperbolicCoolingTower(H_inlet, D_outlet, H_outlet, D_inlet=None, D_base=None, D_throat=None, H_throat=None, H_support=None, D_support=None, n_support=None, inlet_rounding=None)[source]

Bases: object

Class representing the geometry of a hyperbolic cooling tower, as used in many industries especially the poewr industry. All parameters are also attributes.

H_inlet, D_outlet, and H_outlet are always required. Additionally, one set of the following parameters is required; H_support, D_support, n_support, and inlet_rounding are all optional as well.

  • Inlet diameter

  • Inlet diameter and throat diameter

  • Inlet diameter and throat height

  • Inlet diameter, throat diameter, and throat height

  • Base diameter, throat diameter, and throat height

If the inlet diameter is provided but the throat diameter and/or the throat height are missing, two heuristics are used to estimate them (to avoid these heuristics simply specify the values):

  • Assume the throat elevation is 2/3 the elevation of the tower.

  • Assume the throat diameter is 63% the diameter of the inlet.

Parameters:
H_inletpython:float

Height of the inlet zone of the cooling tower (also called rain zone), [m]

D_outletpython:float

The inside diameter of the cooling tower outlet (top of the tower; the elevation the concrete section ends), [m]

H_outletpython:float

The height of the cooling tower outlet (top of the tower;the elevation the concrete section ends), [m]

D_inletpython:float, optional

The inside diameter of the cooling tower inlet at the elevation the concrete section begins, [m]

D_basepython:float, optional

The diameter of the cooling tower at the very base of the tower (the bottom of the inlet zone, at the elevation of the ground), [m]

D_throatpython:float, optional

The diameter of the cooling tower at its minimum section, called its throat; where the two hyperbolas meet, [m]

h_throatpython:float, optional

The elevation of the cooling tower’s throat (its minimum section; where the two hyperbolas meet), [m]

inlet_roundingpython:float, optional

Radius of an optional rounded protrusion from the lip of the cooling tower shell base, which curves upwards from the lip (used to reduce the dead zone area rather than having a flat lip), [m]

H_supportpython:float, optional

The height of each support column, [m]

D_supportpython:float, optional

The diameter of each support column, [m]

n_supportpython:int, optional

The number of support columns of the cooling tower, [m]

Notes

Note there are two hyperbolas in a hyperbolic cooling tower - one under the throat and one above it; they are not necessarily the same.

A hyperbolic cooling tower is not the absolute optimal design, but is is close. The optimality is determined by the amount of material required to build it while maintaining its rigidity. For thermal design purposes, a hyperbolic model covers any minor variation quite well.

References

[1]

Chen, W. F., and E. M. Lui, eds. Handbook of Structural Engineering, Second Edition. Boca Raton, Fla: CRC Press, 2005.

[2]

Ansary, A. M. El, A. A. El Damatty, and A. O. Nassef. Optimum Shape and Design of Cooling Towers, 2011.

Examples

>>> ct = HyperbolicCoolingTower(D_outlet=89.0, H_outlet=200, D_inlet=136.18, H_inlet=14.5)
>>> ct
<Hyperbolic cooling tower, inlet diameter=136.18 m, outlet diameter=89 m, inlet height=14.5 m, outlet height=200 m, throat diameter=85.7934 m, throat height=133.333 m, base diameter=146.427 m>
>>> ct.diameter(5)
142.84514486126062
Attributes:
b_lowerpython:float

The b parameter in the hyperbolic equation for the lower section of the cooling tower, [m]

b_upperpython:float

The b parameter in the hyperbolic equation for the upper section of the cooling tower, [m]

Methods

diameter(H)

Calculates cooling tower diameter at a specified height, using the formulas for either hyperbola, depending on the height specified.

plot

diameter(H)[source]

Calculates cooling tower diameter at a specified height, using the formulas for either hyperbola, depending on the height specified.

\[D = D_{throat}\frac{\sqrt{H^2 + b^2}}{b}\]

The value of H and b used in the above equation is as follows:

  • H_throat - H and b_lower if under the throat

  • H - H_throat and b_upper, if above the throat

Parameters:
Hpython:float

Height at which to calculate the cooling tower diameter, [m]

Returns:
Dpython:float

Diameter of the cooling tower at the specified height, [m]

plot(pts=100)[source]
class fluids.geometry.RectangularFinExchanger(fin_height, fin_thickness, fin_spacing, length=None, width=None, layers=None, plate_thickness=None, flow='crossflow')[source]

Bases: object

Class representing a plate-fin heat exchanger with straight rectangular fins. All parameters are also attributes.

Parameters:
fin_heightpython:float

The total distance between the two metal plates sandwiching the fins and holding them together (abbreviated h), [m]

fin_thicknesspython:float

The thickness of the material the fins were formed from (abbreviated t), [m]

fin_spacingpython:float

The unit cell spacing from one fin to the next; the space between the sides of two fins plus one thickness (abbreviated s), [m]

lengthpython:float, optional

The total length of the flow passage of the plate-fin exchanger (abbreviated L), [m]

widthpython:float, optional

The total width of the space the fins are in; this is also \(N_{fins}\times s\) (abbreviated W), [m]

layerspython:int, optional

The number of layers in the plate-fin exchanger; note these HX almost always single-pass only, [-]

plate_thicknesspython:float, optional

The thickness of the metal separator between layers, [m]

flowpython:str, optional

One of ‘counterflow’, ‘crossflow’, or ‘parallelflow’

Notes

The only required parameters are the fin geometry itself; fin_height, fin_thickness, and fin_spacing.

References

[1]

Yang, Yujie, and Yanzhong Li. “General Prediction of the Thermal Hydraulic Performance for Plate-Fin Heat Exchanger with Offset Strip Fins.” International Journal of Heat and Mass Transfer 78 (November 1, 2014): 860-70. doi:10.1016/j.ijheatmasstransfer.2014.07.060.

[2]

Sheik Ismail, L., R. Velraj, and C. Ranganayakulu. “Studies on Pumping Power in Terms of Pressure Drop and Heat Transfer Characteristics of Compact Plate-Fin Heat Exchangers-A Review.” Renewable and Sustainable Energy Reviews 14, no. 1 (January 2010): 478-85. doi:10.1016/j.rser.2009.06.033.

Examples

>>> PFE = RectangularFinExchanger(0.03, 0.001, 0.012)
>>> PFE.Dh
0.01595
Attributes:
channel_heightpython:float

The height of the channel the fluid flows in \(\text{channel height } = \text{fin height} - \text{fin thickness}\), [m]

channel_widthpython:float

The width of the channel the fluid flows in \(\text{channel width } = \text{fin spacing} - \text{fin thickness}\), [m]

fin_countpython:int

The number of fins per unit length of the layer, \(\text{fin count} = \frac{1}{\text{fin spacing}}\), [1/m]

blockage_ratiopython:float

The fraction of the layer which is blocked to flow by the fins, \(\text{blockage ratio} = \frac{s\cdot h - s\cdot t - t(h-t)}{s\cdot h}\), [m]

A_channelpython:float

Flow area of a single channel in a single layer, \(\text{channel area} = (s-t)(h-t)\), [m]

P_channelpython:float

Wetted perimeter of a single channel in a single layer, \(\text{channel perimeter} = 2(s-t) + 2(h-t)\), [m]

Dhpython:float

Hydraulic diameter of a single channel in a single layer, \(D_{hydraulic} = \frac{4 A_{channel}}{P_{channel}}\), [m]

layer_thicknesspython:float

The thickness of a single layer - the sum of a fin height and a plate thickness, [m]

layer_fin_countpython:int

The number of fins in a layer; rounded to the nearest whole fin, [-]

A_HX_layerpython:float

The surface area including fins for heat transfer in one layer of the HX, [m^2]

A_HXpython:float

The total surface area of the heat exchanger with all layers combined, [m^2]

heightpython:float

The height of all the layers of the heat exchanger combined, plus one extra plate thickness, [m]

volumepython:float

The product of the height, width, and length of the HX, [m^3]

A_specific_HXpython:float

The specific surface area of the heat exchanger - square meters per meter cubed, [m^3]

Methods

set_overall_geometry

set_overall_geometry()[source]
class fluids.geometry.RectangularOffsetStripFinExchanger(fin_length, fin_height, fin_thickness, fin_spacing, length=None, width=None, layers=None, plate_thickness=None, flow='crossflow')[source]

Bases: RectangularFinExchanger

Methods

set_overall_geometry

Tank Volume Functions

fluids.geometry.V_partial_sphere(D, h)[source]

Calculates volume of a partial sphere according to [1]. If h is half of D, the shape is half a sphere. No bottom is considered in this function. Valid inputs are positive values of D and h, with h always smaller or equal to D.

\[a = \sqrt{h(2r - h)}\]
\[V = 1/6 \pi h(3a^2 + h^2)\]
Parameters:
Dpython:float

Diameter of the sphere, [m]

hpython:float

Height, as measured up to where the sphere is cut off, [m]

Returns:
Vpython:float

Volume [m^3]

References

[1]

Weisstein, Eric W. “Spherical Cap.” Text. Accessed December 22, 2015. http://mathworld.wolfram.com/SphericalCap.html.

Examples

>>> V_partial_sphere(1., 0.7)
0.4105014400690663
fluids.geometry.V_horiz_conical(D, L, a, h, headonly=False)[source]

Calculates volume of a tank with conical ends, according to [1].

\[\begin{split}V_f = A_fL + \frac{2aR^2}{3}K, \;\;0 \le h < R\\\end{split}\]
\[\begin{split}V_f = A_fL + \frac{2aR^2}{3}\pi/2,\;\; h = R\\\end{split}\]
\[V_f = A_fL + \frac{2aR^2}{3}(\pi-K), \;\; R< h \le 2R\]
\[K = \cos^{-1} M + M^3\cosh^{-1} \frac{1}{M} - 2M\sqrt{1 - M^2}\]
\[M = \left|\frac{R-h}{R}\right|\]
\[A_f = R^2\cos^{-1}\frac{R-h}{R} - (R-h)\sqrt{2Rh - h^2}\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

Lpython:float

Length of the main cylindrical section, [m]

apython:float

Distance the cone head extends on one side, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

headonlybool, optional

Function returns only the volume of a single head side if True

Returns:
Vpython:float

Volume [m^3]

References

[1] (1,2)

Jones, D. “Calculating Tank Volume.” Text. Accessed December 22, 2015. http://www.webcalc.com.br/blog/Tank_Volume.PDF

Examples

Matching example from [1], with inputs in inches and volume in gallons.

>>> V_horiz_conical(D=108., L=156., a=42., h=36)/231
2041.1923581273443
fluids.geometry.V_horiz_ellipsoidal(D, L, a, h, headonly=False)[source]

Calculates volume of a tank with ellipsoidal ends, according to [1].

\[V_f = A_fL + \pi a h^2\left(1 - \frac{h}{3R}\right)\]
\[A_f = R^2\cos^{-1}\frac{R-h}{R} - (R-h)\sqrt{2Rh - h^2}\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

Lpython:float

Length of the main cylindrical section, [m]

apython:float

Distance the ellipsoidal head extends on one side, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

headonlybool, optional

Function returns only the volume of a single head side if True

Returns:
Vpython:float

Volume [m^3]

References

[1] (1,2)

Jones, D. “Calculating Tank Volume.” Text. Accessed December 22, 2015. http://www.webcalc.com.br/blog/Tank_Volume.PDF

Examples

Matching example from [1], with inputs in inches and volume in gallons.

>>> V_horiz_ellipsoidal(D=108, L=156, a=42, h=36)/231.
2380.9565415578145
fluids.geometry.V_horiz_guppy(D, L, a, h, headonly=False)[source]

Calculates volume of a tank with guppy heads, according to [1].

\[V_f = A_fL + \frac{2aR^2}{3}\cos^{-1}\left(1 - \frac{h}{R}\right) +\frac{2a}{9R}\sqrt{2Rh - h^2}(2h-3R)(h+R)\]
\[A_f = R^2\cos^{-1}\frac{R-h}{R} - (R-h)\sqrt{2Rh - h^2}\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

Lpython:float

Length of the main cylindrical section, [m]

apython:float

Distance the guppy head extends on one side, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

headonlybool, optional

Function returns only the volume of a single head side if True

Returns:
Vpython:float

Volume [m^3]

References

[1] (1,2)

Jones, D. “Calculating Tank Volume.” Text. Accessed December 22, 2015. http://www.webcalc.com.br/blog/Tank_Volume.PDF

Examples

Matching example from [1], with inputs in inches and volume in gallons.

>>> V_horiz_guppy(D=108., L=156., a=42., h=36)/231.
1931.7208029476762
fluids.geometry.V_horiz_spherical(D, L, a, h, headonly=False)[source]

Calculates volume of a tank with spherical heads, according to [1].

\[V_f = A_fL + \frac{\pi a}{6}(3R^2 + a^2),\;\; h = R, |a|\le R\]
\[V_f = A_fL + \frac{\pi a}{3}(3R^2 + a^2),\;\; h = D, |a|\le R\]
\[V_f = A_fL + \pi a h^2\left(1 - \frac{h}{3R}\right),\;\; h = 0, \text{ or } |a| = 0, R, -R\]
\[V_f = A_fL + \frac{a}{|a|}\left\{\frac{2r^3}{3}\left[\cos^{-1} \frac{R^2 - rw}{R(w-r)} + \cos^{-1}\frac{R^2 + rw}{R(w+r)} - \frac{z}{r}\left(2 + \left(\frac{R}{r}\right)^2\right) \cos^{-1}\frac{w}{R}\right] - 2\left(wr^2 - \frac{w^3}{3}\right) \tan^{-1}\frac{y}{z} + \frac{4wyz}{3}\right\} ,\;\; h \ne R, D; a \ne 0, R, -R, |a| \ge 0.01D\]
\[V_f = A_fL + \frac{a}{|a|}\left[2\int_w^R(r^2 - x^2)\tan^{-1} \sqrt{\frac{R^2-x^2}{r^2-R^2}}dx - A_f z\right] ,\;\; h \ne R, D; a \ne 0, R, -R, |a| < 0.01D\]
\[A_f = R^2\cos^{-1}\frac{R-h}{R} - (R-h)\sqrt{2Rh - h^2}\]
\[r = \frac{a^2 + R^2}{2|a|}\]
\[w = R - h\]
\[y = \sqrt{2Rh-h^2}\]
\[z = \sqrt{r^2 - R^2}\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

Lpython:float

Length of the main cylindrical section, [m]

apython:float

Distance the spherical head extends on one side, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

headonlybool, optional

Function returns only the volume of a single head side if True

Returns:
Vpython:float

Volume [m^3]

References

[1] (1,2)

Jones, D. “Calculating Tank Volume.” Text. Accessed December 22, 2015. http://www.webcalc.com.br/blog/Tank_Volume.PDF

Examples

Matching example from [1], with inputs in inches and volume in gallons.

>>> V_horiz_spherical(D=108., L=156., a=42., h=36)/231.
2303.96151169861
fluids.geometry.V_horiz_torispherical(D, L, f, k, h, headonly=False)[source]

Calculates volume of a tank with torispherical heads, according to [1].

\[\begin{split}V_f = A_fL + 2V_1, \;\; 0 \le h \le h_1\\ V_f = A_fL + 2(V_{1,max} + V_2 + V_3), \;\; h_1 < h < h_2\\ V_f = A_fL + 2[2V_{1,max} - V_1(h=D-h) + V_{2,max} + V_{3,max}] , \;\; h_2 \le h \le D\end{split}\]
\[V_1 = \int_0^{\sqrt{2kDh - h^2}} \left[n^2\sin^{-1}\frac{\sqrt {n^2-w^2}}{n} - w\sqrt{n^2-w^2}\right]dx\]
\[V_2 = \int_0^{kD\cos\alpha}\left[n^2\left(\cos^{-1}\frac{w}{n} - \cos^{-1}\frac{g}{n}\right) - w\sqrt{n^2 - w^2} + g\sqrt{n^2 - g^2}\right]dx\]
\[V_3 = \int_w^g(r^2 - x^2)\tan^{-1}\frac{\sqrt{g^2 - x^2}}{z}dx - \frac{z}{2}\left(g^2\cos^{-1}\frac{w}{g} - w\sqrt{2g(h-h_1) - (h-h_1)^2}\right)\]
\[V_{1,max} = v_1(h=h_1)\]
\[v_{2,max} = v_2(h=h_2)\]
\[v_{3,max} = \frac{\pi a_1}{6}(3g^2 + a_1^2)\]
\[a_1 = fD(1-\cos\alpha)\]
\[\alpha = \sin^{-1}\frac{1-2k}{2(f-k)}\]
\[n = R - kD + \sqrt{k^2D^2-x^2}\]
\[g = r\sin\alpha\]
\[r = fD\]
\[h_2 = D - h_1\]
\[w = R - h\]
\[z = \sqrt{r^2- g^2}\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

Lpython:float

Length of the main cylindrical section, [m]

fpython:float

Dimensionless dish-radius parameter; also commonly given as the product of f and D (fD), which is called both dish radius and also crown radius and has units of length, [-]

kpython:float

Dimensionless knuckle-radius parameter; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

hpython:float

Height, as measured up to where the fluid ends, [m]

headonlybool, optional

Function returns only the volume of a single head side if True

Returns:
Vpython:float

Volume [m^3]

References

[1] (1,2)

Jones, D. “Calculating Tank Volume.” Text. Accessed December 22, 2015. http://www.webcalc.com.br/blog/Tank_Volume.PDF

Examples

Matching example from [1], with inputs in inches and volume in gallons.

>>> V_horiz_torispherical(D=108., L=156., f=1., k=0.06, h=36)/231.
2028.62667
fluids.geometry.V_vertical_conical(D, a, h)[source]

Calculates volume of a vertical tank with a convex conical bottom, according to [1]. No provision for the top of the tank is made here.

\[V_f = \frac{\pi}{4}\left(\frac{Dh}{a}\right)^2\left(\frac{h}{3}\right),\; h < a\]
\[V_f = \frac{\pi D^2}{4}\left(h - \frac{2a}{3}\right),\; h\ge a\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Distance the cone head extends under the main cylinder, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

Returns:
Vpython:float

Volume [m^3]

References

[1] (1,2)

Jones, D. “Calculating Tank Volume.” Text. Accessed December 22, 2015. http://www.webcalc.com.br/blog/Tank_Volume.PDF

Examples

Matching example from [1], with inputs in inches and volume in gallons.

>>> V_vertical_conical(132., 33., 24)/231.
250.67461381371024
fluids.geometry.V_vertical_ellipsoidal(D, a, h)[source]

Calculates volume of a vertical tank with a convex ellipsoidal bottom, according to [1]. No provision for the top of the tank is made here.

\[V_f = \frac{\pi}{4}\left(\frac{Dh}{a}\right)^2 \left(a - \frac{h}{3}\right),\; h < a\]
\[V_f = \frac{\pi D^2}{4}\left(h - \frac{a}{3}\right),\; h \ge a\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Distance the ellipsoid head extends under the main cylinder, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

Returns:
Vpython:float

Volume [m^3]

References

[1] (1,2)

Jones, D. “Calculating Tank Volume.” Text. Accessed December 22, 2015. http://www.webcalc.com.br/blog/Tank_Volume.PDF

Examples

Matching example from [1], with inputs in inches and volume in gallons.

>>> V_vertical_ellipsoidal(132., 33., 24)/231.
783.3581681678445
fluids.geometry.V_vertical_spherical(D, a, h)[source]

Calculates volume of a vertical tank with a convex spherical bottom, according to [1]. No provision for the top of the tank is made here.

\[V_f = \frac{\pi h^2}{4}\left(2a + \frac{D^2}{2a} - \frac{4h}{3}\right),\; h < a\]
\[V_f = \frac{\pi}{4}\left(\frac{2a^3}{3} - \frac{aD^2}{2} + hD^2\right),\; h\ge a\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Distance the spherical head extends under the main cylinder, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

Returns:
Vpython:float

Volume [m^3]

References

[1] (1,2)

Jones, D. “Calculating Tank Volume.” Text. Accessed December 22, 2015. http://www.webcalc.com.br/blog/Tank_Volume.PDF

Examples

Matching example from [1], with inputs in inches and volume in gallons.

>>> V_vertical_spherical(132., 33., 24)/231.
583.6018352850442
fluids.geometry.V_vertical_torispherical(D, f, k, h)[source]

Calculates volume of a vertical tank with a convex torispherical bottom, according to [1]. No provision for the top of the tank is made here.

\[V_f = \frac{\pi h^2}{4}\left(2a_1 + \frac{D_1^2}{2a_1} - \frac{4h}{3}\right),\; 0 \le h \le a_1\]
\[V_f = \frac{\pi}{4}\left(\frac{2a_1^3}{3} + \frac{a_1D_1^2}{2}\right) +\pi u\left[\left(\frac{D}{2}-kD\right)^2 +s\right] + \frac{\pi tu^2}{2} - \frac{\pi u^3}{3} + \pi D(1-2k)\left[ \frac{2u-t}{4}\sqrt{s+tu-u^2} + \frac{t\sqrt{s}}{4} + \frac{k^2D^2}{2}\left(\cos^{-1}\frac{t-2u}{2kD}-\alpha\right)\right] ,\; a_1 < h \le a_1 + a_2\]
\[V_f = \frac{\pi}{4}\left(\frac{2a_1^3}{3} + \frac{a_1D_1^2}{2}\right) +\frac{\pi t}{2}\left[\left(\frac{D}{2}-kD\right)^2 +s\right] +\frac{\pi t^3}{12} + \pi D(1-2k)\left[\frac{t\sqrt{s}}{4} + \frac{k^2D^2}{2}\sin^{-1}(\cos\alpha)\right] + \frac{\pi D^2}{4}[h-(a_1+a_2)] ,\; a_1 + a_2 < h\]
\[\alpha = \sin^{-1}\frac{1-2k}{2(f-k)}\]
\[a_1 = fD(1-\cos\alpha)\]
\[a_2 = kD\cos\alpha\]
\[D_1 = 2fD\sin\alpha\]
\[s = (kD\sin\alpha)^2\]
\[t = 2a_2\]
\[u = h - fD(1-\cos\alpha)\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

fpython:float

Dimensionless dish-radius parameter; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

kpython:float

Dimensionless knuckle-radius parameter; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

hpython:float

Height, as measured up to where the fluid ends, [m]

Returns:
Vpython:float

Volume [m^3]

References

[1] (1,2)

Jones, D. “Calculating Tank Volume.” Text. Accessed December 22, 2015. http://www.webcalc.com.br/blog/Tank_Volume.PDF

Examples

Matching example from [1], with inputs in inches and volume in gallons.

>>> V_vertical_torispherical(D=132., f=1.0, k=0.06, h=24)/231.
904.0688283793
fluids.geometry.V_vertical_conical_concave(D, a, h)[source]

Calculates volume of a vertical tank with a concave conical bottom, according to [1]. No provision for the top of the tank is made here.

\[V = \frac{\pi D^2}{12} \left(3h + a - \frac{(a+h)^3}{a^2}\right) ,\;\; 0 \le h < |a|\]
\[V = \frac{\pi D^2}{12} (3h + a ),\;\; h \ge |a|\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Negative distance the cone head extends inside the main cylinder, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

Returns:
Vpython:float

Volume [m^3]

References

[1] (1,2)

Jones, D. “Compute Fluid Volumes in Vertical Tanks.” Chemical Processing. December 18, 2003. http://www.chemicalprocessing.com/articles/2003/193/

Examples

Matching example from [1], with inputs in inches and volume in gallons.

>>> V_vertical_conical_concave(D=113., a=-33, h=15)/231
251.158255657951
fluids.geometry.V_vertical_ellipsoidal_concave(D, a, h)[source]

Calculates volume of a vertical tank with a concave ellipsoidal bottom, according to [1]. No provision for the top of the tank is made here.

\[V = \frac{\pi D^2}{12} \left(3h + 2a - \frac{(a+h)^2(2a-h)}{a^2}\right) ,\;\; 0 \le h < |a|\]
\[V = \frac{\pi D^2}{12} (3h + 2a ),\;\; h \ge |a|\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Negative distance the eppilsoid head extends inside the main cylinder, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

Returns:
Vpython:float

Volume [m^3]

References

[1] (1,2)

Jones, D. “Compute Fluid Volumes in Vertical Tanks.” Chemical Processing. December 18, 2003. http://www.chemicalprocessing.com/articles/2003/193/

Examples

Matching example from [1], with inputs in inches and volume in gallons.

>>> V_vertical_ellipsoidal_concave(D=113., a=-33, h=15)/231
44.84968851034856
fluids.geometry.V_vertical_spherical_concave(D, a, h)[source]

Calculates volume of a vertical tank with a concave spherical bottom, according to [1]. No provision for the top of the tank is made here.

\[V = \frac{\pi}{12}\left[3D^2h + \frac{a}{2}(3D^2 + 4a^2) + (a+h)^3 \left(4 - \frac{3D^2 + 12a^2}{2a(a+h)}\right)\right],\;\; 0 \le h < |a|\]
\[V = \frac{\pi}{12}\left[3D^2h + \frac{a}{2}(3D^2 + 4a^2) \right] ,\;\; h \ge |a|\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Negative distance the spherical head extends inside the main cylinder, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

Returns:
Vpython:float

Volume [m^3]

References

[1] (1,2)

Jones, D. “Compute Fluid Volumes in Vertical Tanks.” Chemical Processing. December 18, 2003. http://www.chemicalprocessing.com/articles/2003/193/

Examples

Matching example from [1], with inputs in inches and volume in gallons.

>>> V_vertical_spherical_concave(D=113., a=-33, h=15)/231
112.81405437348528
fluids.geometry.V_vertical_torispherical_concave(D, f, k, h)[source]

Calculates volume of a vertical tank with a concave torispherical bottom, according to [1]. No provision for the top of the tank is made here.

\[V = \frac{\pi D^2 h}{4} - v_1(h=a_1+a_2) + v_1(h=a_1 + a_2 -h),\; 0 \le h < a_2\]
\[V = \frac{\pi D^2 h}{4} - v_1(h=a_1+a_2) + v_2(h=a_1 + a_2 -h),\; a_2 \le h < a_1 + a_2\]
\[V = \frac{\pi D^2 h}{4} - v_1(h=a_1+a_2) + 0,\; h \ge a_1 + a_2\]
\[v_1 = \frac{\pi}{4}\left(\frac{2a_1^3}{3} + \frac{a_1D_1^2}{2}\right) +\pi u\left[\left(\frac{D}{2}-kD\right)^2 +s\right] + \frac{\pi tu^2}{2} - \frac{\pi u^3}{3} + \pi D(1-2k)\left[ \frac{2u-t}{4}\sqrt{s+tu-u^2} + \frac{t\sqrt{s}}{4} + \frac{k^2D^2}{2}\left(\cos^{-1}\frac{t-2u}{2kD}-\alpha\right)\right]\]
\[v_2 = \frac{\pi h^2}{4}\left(2a_1 + \frac{D_1^2}{2a_1} - \frac{4h}{3}\right)\]
\[\alpha = \sin^{-1}\frac{1-2k}{2(f-k)}\]
\[a_1 = fD(1-\cos\alpha)\]
\[a_2 = kD\cos\alpha\]
\[D_1 = 2fD\sin\alpha\]
\[s = (kD\sin\alpha)^2\]
\[t = 2a_2\]
\[u = h - fD(1-\cos\alpha)\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

fpython:float

Dimensionless dish-radius parameter; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

kpython:float

Dimensionless knuckle-radius parameter; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

hpython:float

Height, as measured up to where the fluid ends, [m]

Returns:
Vpython:float

Volume [m^3]

References

[1] (1,2)

Jones, D. “Compute Fluid Volumes in Vertical Tanks.” Chemical Processing. December 18, 2003. http://www.chemicalprocessing.com/articles/2003/193/

Examples

Matching example from [1], with inputs in inches and volume in gallons.

>>> V_vertical_torispherical_concave(D=113., f=0.71, k=0.081, h=15)/231
103.88569287163769

Tank Surface Area Functions

fluids.geometry.SA_partial_sphere(D, h)[source]

Calculates surface area of a partial sphere according to [1]. If h is half of D, the shape is half a sphere. No bottom is considered in this function. Valid inputs are positive values of D and h, with h always smaller or equal to D.

\[a = \sqrt{h(2r - h)}\]
\[A = \pi(a^2 + h^2)\]
Parameters:
Dpython:float

Diameter of the sphere, [m]

hpython:float

Height, as measured from the cap to where the sphere is cut off [m]

Returns:
SApython:float

Surface area [m^2]

References

[1]

Weisstein, Eric W. “Spherical Cap.” Text. Accessed December 22, 2015. http://mathworld.wolfram.com/SphericalCap.html.

Examples

>>> SA_partial_sphere(1., 0.7)
2.199114857512855
fluids.geometry.SA_ellipsoidal_head(D, a)[source]

Calculates the surface area of an ellipsoidal head according to [1] and [2]. The formula below is for the full shape, the result of which is halved. The formula is for \(a < R\). In the equations, a is the same and c is D.

\[\text{SA} = 2\pi a^2 + \frac{\pi c^2}{e_1}\ln\left(\frac{1+e_1}{1-e_1} \right)\]
\[e_1 = \sqrt{1 - \frac{c^2}{a^2}}\]

For the case of \(a \ge R\) from [2], which is needed to make the tank head volume grow linearly with length:

\[\text{SA} = 2\pi R^2 + \frac{2\pi a^2 R}{\sqrt{a^2 - R^2}}\cos^{-1}\frac{R}{|a|}\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Distance the ellipsoidal head extends, [m]

Returns:
SApython:float

Surface area [m^2]

References

[1]

Weisstein, Eric W. “Spheroid.” Text. Accessed March 14, 2016. http://mathworld.wolfram.com/Spheroid.html.

[2] (1,2)

Jones, D. “Calculating Tank Wetted Area.” Text. Chemical Processing. April 2017. https://www.chemicalprocessing.com/assets/Uploads/calculating-tank-wetted-area.pdf

Examples

Spherical case

>>> SA_ellipsoidal_head(2, 1)
6.283185307179586
>>> SA_ellipsoidal_head(2, 1.5)
8.459109081729984
fluids.geometry.SA_conical_head(D, a)[source]

Calculates the surface area of a conical head according to [1].

\[SA = \frac{\pi D}{2} \sqrt{a^2 + \left(\frac{D}{2}\right)^2}\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Distance the conical head extends, [m]

Returns:
SApython:float

Surface area [m^2]

References

[1]

Weisstein, Eric W. “Cone.” Text. Accessed March 14, 2016. http://mathworld.wolfram.com/Cone.html.

Examples

>>> SA_conical_head(2, 1)
4.442882938158366
fluids.geometry.SA_guppy_head(D, a)[source]

Calculates the surface area of a guppy head according to [1]. Some work was involved in combining formulas for the ellipse of the head, and the conic section on the sides.

\[SA = \frac{\pi D}{4}\sqrt{D^2 + a^2} + \frac{\pi D}{2}a\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Distance the conical head extends, [m]

Returns:
SApython:float

Surface area [m^2]

References

[1]

Weisstein, Eric W. “Cone.” Text. Accessed March 14, 2016. http://mathworld.wolfram.com/Cone.html.

Examples

>>> SA_guppy_head(2, 1)
6.654000019110157
fluids.geometry.SA_torispheroidal(D, f, k)[source]

Calculates surface area of a torispherical head according to [1]. Somewhat involved. Equations are adapted to be used for a full head.

\[SA = S_1 + S_2\]
\[S_1 = 2\pi D^2 f_d \alpha\]
\[S_2 = 2\pi D^2 f_k\left(\alpha - \alpha_1 + (0.5 - f_k)\left(\sin^{-1} \left(\frac{\alpha-\alpha_2}{f_k}\right) - \sin^{-1}\left(\frac{ \alpha_1-\alpha_2}{f_k}\right)\right)\right)\]
\[\alpha_1 = f_d\left(1 - \sqrt{1 - \left(\frac{0.5 - f_k}{f_d-f_k} \right)^2}\right)\]
\[\alpha_2 = f_d - \sqrt{f_d^2 - 2f_d f_k + f_k - 0.25}\]
\[\alpha = \frac{a}{D_i}\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

fpython:float

Dimensionless dish-radius parameter; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

kpython:float

Dimensionless knuckle-radius parameter; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

Returns:
SApython:float

Surface area [m^2]

References

[1] (1,2)

Honeywell. “Calculate Surface Areas and Cross-sectional Areas in Vessels with Dished Heads”. https://www.honeywellprocess.com/library/marketing/whitepapers/WP-VesselsWithDishedHeads-UniSimDesign.pdf Whitepaper. 2014.

Examples

Example from [1].

>>> SA_torispheroidal(D=2.54, f=1.039370079, k=0.062362205)
6.00394283477063
fluids.geometry.SA_partial_cylindrical_body(L, D, h)[source]

Calculates the partial area of a cylinder’s body in the context of a horizontal cylindrical vessel and liquid partially filling it. This computes the wetted surface area of the bottom of the cylinder.

\[\text{SA} = L D \cos^{-1}\left(\frac{D - 2h}{D}\right)\]
Parameters:
Lpython:float

Length of the cylinder, [m]

Dpython:float

Diameter of the cylinder, [m]

hpython:float

Height measured from bottom of cylinder to liquid level, [m]

Returns:
SA_partialpython:float

Partial (wetted) surface area, [m^2]

Notes

This method is undefined for \(h > D\). and \(h < 0\), but those cases are handled by returning the full surface area and the zero respectively.

References

[1]

Weisstein, Eric W. “Circular Segment.” Text. Wolfram Research, Inc. Accessed May 10, 2020. https://mathworld.wolfram.com/CircularSegment.html.

Examples

>>> SA_partial_cylindrical_body(L=200.0, D=96., h=22.0)
19168.852890279868
fluids.geometry.SA_partial_horiz_conical_head(D, a, h)[source]

Calculates the partial area of a conical tank head in the context of a horizontal vessel and liquid partially filling it. This computes the wetted surface area of one of the conical heads only.

\[\text{SA} = \frac{\sqrt{(a^2 + R^2)}}{R}\left[R^2\cos^{-1}\left(\frac{ (R-h)}{R}\right) - (R-h)\sqrt{(2Rh - h^2)}\right]\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Distance the cone head extends on one side, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

Returns:
SA_partialpython:float

Partial (wetted) surface area of one conical tank head, [m^2]

Notes

This method is undefined for \(h > D\) and \(h < 0\), but those cases are handled by returning the full surface area and the zero respectively.

References

[1]

Jones, D. “Calculating Tank Wetted Area.” Text. Chemical Processing. April 2017. https://www.chemicalprocessing.com/assets/Uploads/calculating-tank-wetted-area.pdf

Examples

>>> SA_partial_horiz_conical_head(D=72., a=48.0, h=24.0)
1980.0498315169873
fluids.geometry.SA_partial_horiz_spherical_head(D, a, h)[source]

Calculates the partial area of a spherical tank head in the context of a horizontal vessel and liquid partially filling it. This computes the wetted surface area of one of the spherical heads only.

\[\text{SA} = \frac{a^2 + R^2}{|a|}\int_{R-h}^R \sin^{-1} \frac{2|a|\sqrt{R^2-x^2}} {\sqrt{(a^2+R^2)^2 - (2ax)^2}} dx\]

For the special case of \(|a| = R\) :

\[\text{SA} = \pi R h\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Distance the spherical head extends on one side, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

Returns:
SA_partialpython:float

Partial (wetted) surface area of one spherical tank head, [m^2]

Notes

This method is undefined for \(h > D\) and \(h < 0\), but those cases are handled by returning the full surface area and the zero respectively.

A symbolic attempt did not suggest any analytical integrals were available.

References

[1]

Jones, D. “Calculating Tank Wetted Area.” Text. Chemical Processing. April 2017. https://www.chemicalprocessing.com/assets/Uploads/calculating-tank-wetted-area.pdf

Examples

>>> SA_partial_horiz_spherical_head(D=72., a=48.0, h=24.0)
2027.2672
fluids.geometry.SA_partial_horiz_ellipsoidal_head(D, a, h)[source]

Calculates the partial area of a ellipsoidal tank head in the context of a horizontal vessel and liquid partially filling it. This computes the wetted surface area of one of the ellipsoidal heads only.

\[\text{SA} = \frac{2}{R} \int_{R-h}^R \int_0^{\sqrt{R^2 - x^2}} \sqrt{ \frac{(R^2 - a^2)x^2 + (R^2 - a^2)y^2 - R^4} {x^2 + y^2 - R^2}} dy dx\]

After extensive manipulation, the first integral was solved analytically, extending the result of [1] with greater performance.

\[\text{SA} = \frac{2}{R} \int_{R-h}^R \frac{\left(\frac{R^{4} - R^{2} \left(R^{2} - a^{2}\right)}{R^{2} - y^{2}}\right)^{0.5} \left(R^{2} - y^{2}\right)^{0.5} E{\left(\frac{\left(- R^{2} + y^{2}\right) \left(R^{2} - a^{2}\right)}{- R^{4} + y^{2} \left(R^{2} - a^{2}\right)} \right)}} {\left(\frac{R^{4} - R^{2} \left(R^{2} - a^{2}\right)}{R^{4} - y^{2} \left(R^{2} - a^{2}\right)}\right)^{0.5}}\]

Where \(E(x)\) is the complete elliptic integral of the second kind, calculated with SciPy’s link to the cephes library.

Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Distance the ellipsoidal head extends on one side, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

Returns:
SA_partialpython:float

Partial (wetted) surface area of one ellipsoidal tank head, [m^2]

Notes

This method is undefined for \(h > D\) and \(h < 0\), but those cases are handled by returning the full surface area and the zero respectively.

The original numerical double integral is extremely nasty - there are places where f(x) -> infinity but that have a bounded area. quadpack’s numerical integration handles this well, but adaptive inetgration which is not aware of singularities does not.

References

[1]

Jones, D. “Calculating Tank Wetted Area.” Text. Chemical Processing. April 2017. https://www.chemicalprocessing.com/assets/Uploads/calculating-tank-wetted-area.pdf

Examples

>>> SA_partial_horiz_ellipsoidal_head(D=72., a=48.0, h=24.0)
3401.233622547
fluids.geometry.SA_partial_horiz_guppy_head(D, a, h)[source]

Calculates the partial area of a guppy tank head in the context of a horizontal vessel and liquid partially filling it. This computes the wetted surface area of one of the guppy heads only.

\[\text{SA} = 2\int_{-R}^{h-R}\int_{0}^{\sqrt{R^2 - x^2}} \sqrt{1 + \left(\frac{a}{2R}\left(1 - \frac{y^2}{(R-x)^2} \right) \right)^2 + \left(\frac{ay}{R(R-x)} \right)^2 } dy dx\]

After extensive manipulation, the first integral was solved analytically, extending the result of [1]. Even with the special functions, this form has somewhat greater performance (and improved precision).

\[\text{SA} = 2 \int_{-R}^{h-R} \frac{\frac{2 a \left(4 + \frac{a^{2} \left(2 R^{2} - 2 R y\right)^{2}}{R^{2} \left(R - y\right)^{4}}\right) \sqrt{R^{2} - y^{2}}}{\sqrt{4 R^{2} + a^{2}} \left(\frac{a \left(R^{2} - y^{2}\right)}{\left(R - y\right)^{2} \sqrt{4 R^{2} + a^{2}}} + 1\right)} + \left(4 + \frac{a^{2} \left(2 R^{2} - 2 R y\right) ^{2}}{R^{2} \left(R - y\right)^{4}}\right) \sqrt{R^{2} - y^{2}} - \frac{2 \sqrt{a} \sqrt{\frac{4 R^{2} \left(R - y\right)^{4} + a^{2} \left(2 R^{2} - 2 R y\right)^{2}}{\left(R^{2} \sqrt{4 R^{2} + a^{2}} - 2 R y \sqrt{4 R^{2} + a^{2}} + a \left(R^{2} - y^{2}\right) + y^{2} \sqrt{4 R^{2} + a^{2}}\right)^{2}}} \left(R - y\right) \left(4 R^{2} + a^{2}\right)^{0.75} \left(\frac{a \left(R^{2} - y^{2}\right)}{\left(R - y\right)^{2} \sqrt{4 R^{2} + a^{2}}} + 1\right) \operatorname{ellipeinc}{\left(2 \operatorname{atan}{\left(\frac{\sqrt{a} \sqrt{R^{2} - y^{2}}} {\left(R - y\right) \left(4 R^{2} + a^{2}\right)^{0.25}} \right)}, - \frac{a}{2 \sqrt{4 R^{2} + a^{2}}} + 0.5 \right)}}{R^{2}} + \frac{1.0 \sqrt{\frac{4 R^{2} \left(R - y\right)^{4} + a^{2} \left(2 R^{2} - 2 R y\right)^{2}}{\left(R^{2} \sqrt{4 R^{2} + a^{2}} - 2 R y \sqrt{4 R^{2} + a^{2}} + a \left(R^{2} - y^{2}\right) + y^{2} \sqrt{4 R^{2} + a^{2}}\right)^{2}}} \left(R - y\right) \left(4 R^{2} + a^{2}\right)^{0.25} \left(\frac{a \left(R^{2} - y^{2}\right)}{\left(R - y\right)^{2} \sqrt{4 R^{2} + a^{2}}} + 1\right) \left(4 R^{2} + a^{2} + a \sqrt{4 R^{2} + a^{2}}\right) \operatorname{ellipkinc}{\left(2 \operatorname{atan}{\left(\frac{\sqrt{a} \sqrt{R^{2} - y^{2}}}{\left(R - y\right) \left(4 R^{2} + a^{2}\right)^{0.25}} \right)},- \frac{a}{2 \sqrt{4 R^{2} + a^{2}}} + 0.5 \right)}}{R^{2} \sqrt{a}}}{6 \sqrt{4 + \frac{a^{2} \left(2 R^{2} - 2 R y\right)^{2}}{R^{2} \left(R - y\right)^{4}}}}\]

Where ellipeinc is the incomplete elliptic integral of the second kind, and ellipkinc is the incomplete elliptic integral of the first kind, both calculated with SciPy’s link to the cephes library.

Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Distance the guppy head extends on one side, [m]

hpython:float

Height, as measured up to where the fluid ends, [m]

Returns:
SA_partialpython:float

Partial (wetted) surface area of one guppy tank head, [m^2]

Notes

This method is undefined for \(h > D\) and \(h < 0\), but those cases are handled by returning the full surface area and the zero respectively.

The analytical integral was derived with Rubi.

References

[1]

Jones, D. “Calculating Tank Wetted Area.” Text. Chemical Processing. April 2017. https://www.chemicalprocessing.com/assets/Uploads/calculating-tank-wetted-area.pdf

Examples

>>> SA_partial_horiz_guppy_head(D=72., a=48.0, h=24.0)
1467.8949
fluids.geometry.SA_partial_horiz_torispherical_head(D, f, k, h)[source]

Calculates the partial area of a torispherical tank head in the context of a horizontal vessel and liquid partially filling it. This computes the wetted surface area of one of the torispherical heads only.

The expressions used are quite complicated; see [1] for more details.

Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

fpython:float

Dimensionless dish-radius parameter; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

kpython:float

Dimensionless knuckle-radius parameter; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

hpython:float

Height, as measured up to where the fluid ends, [m]

Returns:
SA_partialpython:float

Partial (wetted) surface area of one torispherical tank head, [m^2]

Notes

This method is undefined for \(h > D\) and \(h < 0\), but those cases are handled by returning the full surface area and the zero respectively.

One integral:

\[\int_{R-h}^{fD\sin \alpha} \cos^{-1} \frac{fD\cos \alpha}{\sqrt{f^2 D^2 - x^2}} dx\]

Can be computed as follows, using WolframAlpha.

\[x \operatorname{acos}{\left(\frac{c}{\sqrt{b - x^{2}}} \right)} + \frac{\sqrt{b} \sqrt{- \left(b - x^{2}\right)^{2}} \left(- b + c^{2} + x^{2}\right) \operatorname{atan}{\left(\frac{c x}{\sqrt{b} \sqrt{b - c^{2} - x^{2}}} \right)} + c \sqrt{- \left(- b + c^{2} + x^{2}\right)^{2}} \left(- b + x^{2}\right) \operatorname{atan}{\left(\frac{x \sqrt{b - x^{2}}} {\sqrt{- b + x^{2}} \sqrt{- b + c^{2} + x^{2}}} \right)}}{\sqrt{ \frac{- b + c^{2} + x^{2}}{- b + x^{2}}} \left(- b + x^{2}\right)^{1.5} \sqrt{b - c^{2} - x^{2}}}\]

With the following constants:

\[c = fD\cos \alpha\]
\[b = f^2 D^2\]

The other integral is a double integral. There is an analytical integral available for the first integral, which takes the form:

\[2 \sqrt{\frac{R^{2} k^{2} \left(4 R^{2} k^{2} - y^{2} + \left(- 2 R k + R\right)^{2} + 2 \left(- 2 R k + R\right) \sqrt{4 R^{2} k^{2} - y^{2}} \right)}{\left(4 R^{2} k^{2} - y^{2}\right) \left(\left(R - h\right)^{2} - \left(- 4 R k + 2 R\right) \sqrt{4 R^{2} k^{2} - y^{2}} + 2 \left(- 2 R k + R\right) \sqrt{4 R^{2} k^{2} - y^{2}} \right)}} \sqrt{\left(R - h\right)^{2} - \left(- 4 R k + 2 R\right) \sqrt{4 R^{2} k^{2} - y^{2}} + 2 \left(- 2 R k + R\right) \sqrt{4 R^{2} k^{2} - y^{2}}} \operatorname{atan}{\left(\frac{ \sqrt{4 R^{2} k^{2} - y^{2} - \left(R - h\right)^{2} + \left(- 4 R k + 2 R\right) \sqrt{4 R^{2} k^{2} - y^{2}} + \left(- 2 R k + R\right) ^{2}}}{\sqrt{\left(R - h\right)^{2} - \left(- 4 R k + 2 R\right) \sqrt{4 R^{2} k^{2} - y^{2}} + 2 \left(- 2 R k + R\right) \sqrt{4 R^{2} k^{2} - y^{2}}}} \right)}\]

References

[1]

Jones, D. “Calculating Tank Wetted Area.” Text. Chemical Processing. April 2017. https://www.chemicalprocessing.com/assets/Uploads/calculating-tank-wetted-area.pdf

Examples

>>> SA_partial_horiz_torispherical_head(D=72., f=1, k=.06, h=24.0)
1471.201832459
fluids.geometry.SA_partial_vertical_conical_head(D, a, h)[source]

Calculates the partial area of a conical tank head in the context of a vertical vessel and liquid partially filling it. This computes the wetted surface area of one of the conical heads only, and is valid for h up to a only.

\[\text{SA} = \frac{\pi R h^2 \sqrt{a^2 + R^2}}{a^2}\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Distance the cone head extends beneath the vertical tank, [m]

hpython:float

Height, as measured up to where the fluid ends or the top of the conical head, whichever is less, [m]

Returns:
SA_partialpython:float

Partial (wetted) surface area of one conical tank head extending beneath the vessel, [m^2]

Notes

This method is undefined for \(h < 0\), but this is handled by returning zero.

References

[1]

Jones, D. “Calculating Tank Wetted Area.” Text. Chemical Processing. April 2017. https://www.chemicalprocessing.com/assets/Uploads/calculating-tank-wetted-area.pdf

Examples

>>> SA_partial_vertical_conical_head(D=72., a=48.0, h=24.0)
1696.4600329384882
fluids.geometry.SA_partial_vertical_ellipsoidal_head(D, a, h)[source]

Calculates the partial area of a ellipsoidal tank head in the context of a vertical vessel and liquid partially filling it. This computes the wetted surface area of one of the ellipsoidal heads only, and is valid for h up to a only.

If \(a > R\):

\[\text{SA} = \pi R^2 - \frac{\pi (a - h)R\sqrt{a^4 - (a-h)^2(a^2-R^2)}}{a^2} + \frac{\pi a^2 R}{\sqrt{a^2 - R^2}}\left( \cos^{-1} \frac{R}{a} - \sin^{-1} \frac{(a-h)\sqrt{a^2-R^2}}{a^2} \right)\]

Otherwise for \(0 < a < R\):

\[\text{SA} = \pi R^2 - \frac{\pi (a - h)R\sqrt{a^4 - (a-h)^2(a^2-R^2)}}{a^2} + \frac{\pi a^2 R}{\sqrt{a^2 - R^2}}\ln \left(\frac{a(\sqrt{R^2 - a^2} + R)} {(a-h)\sqrt{R^2 - a^2} + \sqrt{a^4 + (a-h)^2(R^2 - a^2)}} \right)\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Distance the ellipsoidal head extends beneath the vertical tank, [m]

hpython:float

Height, as measured up to where the fluid ends or the top of the ellipsoidal head, whichever is less, [m]

Returns:
SA_partialpython:float

Partial (wetted) surface area of one ellipsoidal tank head extending beneath the vessel, [m^2]

Notes

This method is undefined for \(h < 0\), but this is handled by returning zero.

References

[1]

Jones, D. “Calculating Tank Wetted Area.” Text. Chemical Processing. April 2017. https://www.chemicalprocessing.com/assets/Uploads/calculating-tank-wetted-area.pdf

Examples

>>> SA_partial_vertical_ellipsoidal_head(D=72., a=48.0, h=24.0)
4675.23789137632
fluids.geometry.SA_partial_vertical_spherical_head(D, a, h)[source]

Calculates the partial area of a spherical tank head in the context of a vertical vessel and liquid partially filling it. This computes the wetted surface area of one of the conical heads only, and is valid for h up to a only.

\[\text{SA} = \pi h \left(\frac{a^2 + R^2}{a}\right)\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

apython:float

Distance the spherical head extends beneath the vertical tank, [m]

hpython:float

Height, as measured up to where the fluid ends or the top of the spherical head, whichever is less, [m]

Returns:
SA_partialpython:float

Partial (wetted) surface area of one spherical tank head extending beneath the vessel, [m^2]

Notes

This method is undefined for \(h < 0\), but this is handled by returning zero.

References

[1]

Jones, D. “Calculating Tank Wetted Area.” Text. Chemical Processing. April 2017. https://www.chemicalprocessing.com/assets/Uploads/calculating-tank-wetted-area.pdf

Examples

>>> SA_partial_vertical_spherical_head(72, a=24, h=12)
2940.5307237600464
fluids.geometry.SA_partial_vertical_torispherical_head(D, f, k, h)[source]

Calculates the partial area of a torispherical tank head in the context of a vertical vessel and liquid partially filling it. This computes the wetted surface area of one of the torispherical heads only.

if \(a_1 <= h\):

\[\text{SA} = 2\pi f D h\]

if \(a_1 \le h \le a\):

\[\text{SA} = 2\pi f D a_1 + 2\pi k D\left( h - a_1 + (R - kD) \left( \sin^{-1} \frac{a_2}{kD} - \sin^{-1} \frac{a-h}{kD} \right) \right)\]
\[\alpha = \sin^{-1}\frac{1-2k}{2(f-k)}\]
\[a_1 = fD(1-\cos\alpha)\]
\[a_2 = kD\cos\alpha\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

fpython:float

Dimensionless dish-radius parameter; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

kpython:float

Dimensionless knuckle-radius parameter; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

hpython:float

Height, as measured up to where the fluid ends or the top of the torispherical head, whichever is less, [m]

Returns:
SA_partialpython:float

Partial (wetted) surface area of one torispherical tank head, [m^2]

Notes

This method is undefined for \(h > D\) and \(h < 0\), but those cases are handled by returning the full surface area and the zero respectively.

References

[1]

Jones, D. “Calculating Tank Wetted Area.” Text. Chemical Processing. April 2017. https://www.chemicalprocessing.com/assets/Uploads/calculating-tank-wetted-area.pdf

Examples

This method is undefined for \(h < 0\), but this is handled by returning zero.

Miscellaneous Geometry Functions

fluids.geometry.pitch_angle_solver(angle=None, pitch=None, pitch_parallel=None, pitch_normal=None)[source]

Utility to take any two of angle, pitch, pitch_parallel, and pitch_normal and calculate the other two. This is useful for applications with tube banks, as in shell and tube heat exchangers or air coolers and allows for a wider range of user input.

\[\text{pitch normal} = \text{pitch} \cdot \sin(\text{angle})\]
\[\text{pitch parallel} = \text{pitch} \cdot \cos(\text{angle})\]
Parameters:
anglepython:float, optional

The angle of the tube layout, [degrees]

pitchpython:float, optional

The shortest distance between tube centers; defined in relation to the flow direction only, [m]

pitch_parallelpython:float, optional

The distance between tube center along a line parallel to the flow; has been called longitudinal pitch, pp, s2, SL, and p2, [m]

pitch_normalpython:float, optional

The distance between tube centers in a line 90° to the line of flow; has been called the transverse pitch, pn, s1, ST, and p1, [m]

Returns:
anglepython:float

The angle of the tube layout, [degrees]

pitchpython:float

The shortest distance between tube centers; defined in relation to the flow direction only, [m]

pitch_parallelpython:float

The distance between tube center along a line parallel to the flow; has been called longitudinal pitch, pp, s2, SL, and p2, [m]

pitch_normalpython:float

The distance between tube centers in a line 90° to the line of flow; has been called the transverse pitch, pn, s1, ST, and p1, [m]

Notes

For the 90 and 0 degree case, the normal or parallel pitches can be zero; given the angle and the zero value, obviously is it not possible to calculate the pitch and a math error will be raised.

No exception will be raised if three or four inputs are provided; the other two will simply be calculated according to the list of if statements used.

An exception will be raised if only one input is provided.

References

[1]

Schlunder, Ernst U, and International Center for Heat and Mass Transfer. Heat Exchanger Design Handbook. Washington: Hemisphere Pub. Corp., 1983.

Examples

>>> pitch_angle_solver(pitch=1, angle=30)
(30, 1, 0.8660254037844387, 0.49999999999999994)
fluids.geometry.plate_enlargement_factor(amplitude, wavelength)[source]

Calculates the enhancement factor of the sinusoidal waves of the plate heat exchanger. This is the multiplier for the flat plate area to obtain the actual area available for heat transfer. Obtained from the following integral:

\[\phi = \frac{\text{Effective area}}{\text{Projected area}} = \frac{\int_0^\lambda\sqrt{1 + \left(\frac{\gamma\pi}{2}\right)^2 \cos^2\left(\frac{2\pi}{\lambda}x\right)}dx}{\lambda}\]
\[\gamma = \frac{4a}{\lambda}\]

The solution to the integral is:

\[\phi = \frac{2E\left(\frac{-4a^2\pi^2}{\lambda^2}\right)}{\pi}\]

where E is the complete elliptic integral of the second kind, calculated with SciPy.

Parameters:
amplitudepython:float

Half the height of the wave of the ridges, [m]

wavelengthpython:float

Distance between the bottoms of two of the ridges (sometimes called pitch), [m]

Returns:
plate_enlargement_factorpython:float

The extra surface area multiplier as compared to a flat plate caused the corrugations, [-]

Notes

This is the exact analytical integral, obtained via Mathematica, Maple, and quite a bit of trial and error. It is confirmed via numerical integration. The expression normally given is an approximation as follows:

\[\phi = \frac{1}{6}\left(1+\sqrt{1+A^2} + 4\sqrt{1+A^2/2}\right)\]
\[A = \frac{2\pi a}{\lambda}\]

Most plate heat exchangers approximate a sinusoidal geometry only.

Examples

>>> plate_enlargement_factor(amplitude=5E-4, wavelength=3.7E-3)
1.1611862034509677
fluids.geometry.a_torispherical(D, f, k)[source]

Calculates depth of a torispherical head according to [1].

\[a = a_1 + a_2\]
\[\alpha = \sin^{-1}\frac{1-2k}{2(f-k)}\]
\[a_1 = fD(1-\cos\alpha)\]
\[a_2 = kD\cos\alpha\]
Parameters:
Dpython:float

Diameter of the main cylindrical section, [m]

fpython:float

Dimensionless dish-radius parameter; also commonly given as the product of f and D (fD), which is called dish radius and has units of length, [-]

kpython:float

Dimensionless knuckle-radius parameter; also commonly given as the product of k and D (kD), which is called the knuckle radius and has units of length, [-]

Returns:
apython:float

Depth of head [m]

References

[1] (1,2)

Jones, D. “Calculating Tank Volume.” Text. Accessed December 22, 2015. http://www.webcalc.com.br/blog/Tank_Volume.PDF

Examples

Example from [1].

>>> a_torispherical(D=96., f=0.9, k=0.2)
25.684268924767125
fluids.geometry.A_partial_circle(D, h)[source]

Calculates the partial area of a circle, in the context of the circle being an end cap to a horizontal cylindrical vessel and liquid partially filling it. This computes the wetted surface area of one of the end caps.

Multiply this by two to obtain the wetted area of two end caps.

\[\text{SA} = R^2\cos^{-1}\frac{(R - h)}{R} - (R - h)\sqrt{(2Rh - h^2)}\]
Parameters:
Dpython:float

Diameter of the circle, [m]

hpython:float

Height measured from bottom of circle to liquid level, [m]

Returns:
SA_partialpython:float

Partial (wetted) surface area, [m^2]

Notes

This method is undefined for \(h > D\) and \(h < 0\), but those cases are handled by returning the full surface area and the zero respectively.

References

[1]

Weisstein, Eric W. “Circular Segment.” Text. Wolfram Research, Inc. Accessed May 10, 2020. https://mathworld.wolfram.com/CircularSegment.html.

Examples

>>> A_partial_circle(D=96., h=22.0)
1251.2018147383194
fluids.geometry.circle_segment_h_from_A(A, D)[source]

Calculates the height of a chord of a circle given the area of that circle segment. This is a numerical problem, solving the following equation for h.

\[\text{A} = R^2\cos^{-1}\frac{(R - h)}{R} - (R - h)\sqrt{(2Rh - h^2)}\]
Parameters:
Apython:float

Circle section area, [m^2]

Dpython:float

Diameter of the circle, [m]

Returns:
hpython:float

Height measured from bottom of circle to the end of the circle section, [m]

References

[1]

Weisstein, Eric W. “Circular Segment.” Text. Wolfram Research, Inc. Accessed May 10, 2020. https://mathworld.wolfram.com/CircularSegment.html.

Examples

>>> circle_segment_h_from_A(A=1251.2018147383194, D=96.)
22.0

Pellet Properties

fluids.geometry.sphericity(A, V)[source]

Returns the sphericity of a particle of surface area A and volume V. Sphericity is the ratio of the surface area of a sphere with the same volume as the particle (equivalent diameter) to the actual surface area of the particle.

\[\Psi = \frac{\text{A of sphere with } V_p } {{A}_p} = \frac{\pi^{\frac{1}{3}}(6V_p)^{\frac{2}{3}}}{A_p}\]
Parameters:
Apython:float

Surface area of particle, [m^2]

Vpython:float

Volume of particle, [m^3]

Returns:
Psipython:float

Sphericity [-]

Notes

All non-spherical particles have spericities less than 1 but greater than 0. Many common geometrical shapes have their results calculated exactly in [2].

References

[1]

Rhodes, Martin J., ed. Introduction to Particle Technology. 2E. Chichester, England ; Hoboken, NJ: Wiley, 2008.

[2]

“Sphericity.” Wikipedia, March 8, 2017. https://en.wikipedia.org/w/index.php?title=Sphericity&oldid=769183043

Examples

>>> sphericity(10., 2.)
0.767663317071005

For a cube of side length a=3, the surface area is 6*a^2=54 and volume a^3=27. Its sphericity is then:

>>> sphericity(A=54, V=27)
0.8059959770082346
fluids.geometry.aspect_ratio(Dmin, Dmax)[source]

Returns the aspect ratio of a shape with minimum and maximum dimension, Dmin and Dmax.

\[A_R = \frac{D_{min}}{D_{max}}\]
Parameters:
Dminpython:float

Minimum dimension, [m]

Dmaxpython:float

Maximum dimension, [m]

Returns:
a_rpython:float

Aspect ratio [-]

Examples

>>> aspect_ratio(.2, 2)
0.1
fluids.geometry.circularity(A, P)[source]

Returns the circularity of a shape with area A and perimeter P.

\[f_{circ} = \frac {4 \pi A} {P^2}\]

Defined to be 1 for a circle. Used to characterize particles. Any non-circular shape must have a circularity less than one.

Parameters:
Apython:float

Area of the shape, [m^2]

Ppython:float

Perimeter of the shape, [m]

Returns:
f_circpython:float

Circularity of the shape [-]

Examples

Square, side length = 2 (all squares are the same):

>>> circularity(A=(2*2), P=4*2)
0.7853981633974483

Rectangle, one side length = 1, second side length = 100

>>> D1 = 1
>>> D2 = 100
>>> A = D1*D2
>>> P = 2*D1 + 2*D2
>>> circularity(A, P)
0.030796908671598795
fluids.geometry.A_cylinder(D, L)[source]

Returns the surface area of a cylinder.

\[A = \pi D L + 2\cdot \frac{\pi D^2}{4}\]
Parameters:
Dpython:float

Diameter of the cylinder, [m]

Lpython:float

Length of the cylinder, [m]

Returns:
Apython:float

Surface area [m^2]

Examples

>>> A_cylinder(0.01, .1)
0.0032986722862692833
fluids.geometry.V_cylinder(D, L)[source]

Returns the volume of a cylinder.

\[V = \frac{\pi D^2}{4}L\]
Parameters:
Dpython:float

Diameter of the cylinder, [m]

Lpython:float

Length of the cylinder, [m]

Returns:
Vpython:float

Volume [m^3]

Examples

>>> V_cylinder(0.01, .1)
7.853981633974484e-06
fluids.geometry.A_hollow_cylinder(Di, Do, L)[source]

Returns the surface area of a hollow cylinder.

\[A = \pi D_o L + \pi D_i L + 2\cdot \frac{\pi D_o^2}{4} - 2\cdot \frac{\pi D_i^2}{4}\]
Parameters:
Dipython:float

Diameter of the hollow in the cylinder, [m]

Dopython:float

Diameter of the exterior of the cylinder, [m]

Lpython:float

Length of the cylinder, [m]

Returns:
Apython:float

Surface area [m^2]

Examples

>>> A_hollow_cylinder(0.005, 0.01, 0.1)
0.004830198704894308
fluids.geometry.V_hollow_cylinder(Di, Do, L)[source]

Returns the volume of a hollow cylinder.

\[V = \frac{\pi D_o^2}{4}L - L\frac{\pi D_i^2}{4}\]
Parameters:
Dipython:float

Diameter of the hollow in the cylinder, [m]

Dopython:float

Diameter of the exterior of the cylinder, [m]

Lpython:float

Length of the cylinder, [m]

Returns:
Vpython:float

Volume [m^3]

Examples

>>> V_hollow_cylinder(0.005, 0.01, 0.1)
5.890486225480862e-06
fluids.geometry.A_multiple_hole_cylinder(Do, L, holes)[source]

Returns the surface area of a cylinder with multiple holes. Calculation will naively return a negative value or other impossible result if the number of cylinders added is physically impossible. Holes may be of different shapes, but must be perpendicular to the axis of the cylinder.

\[A = \pi D_o L + 2\cdot \frac{\pi D_o^2}{4} + \sum_{i}^n \left( \pi D_i L - 2\cdot \frac{\pi D_i^2}{4}\right)\]
Parameters:
Dopython:float

Diameter of the exterior of the cylinder, [m]

Lpython:float

Length of the cylinder, [m]

holespython:list

List of tuples containing (diameter, count) pairs of descriptions for each of the holes sizes.

Returns:
Apython:float

Surface area [m^2]

Examples

>>> A_multiple_hole_cylinder(0.01, 0.1, [(0.005, 1)])
0.004830198704894308
fluids.geometry.V_multiple_hole_cylinder(Do, L, holes)[source]

Returns the solid volume of a cylinder with multiple cylindrical holes. Calculation will naively return a negative value or other impossible result if the number of cylinders added is physically impossible.

\[V = \frac{\pi D_o^2}{4}L - L\frac{\pi D_i^2}{4}\]
Parameters:
Dopython:float

Diameter of the exterior of the cylinder, [m]

Lpython:float

Length of the cylinder, [m]

holespython:list

List of tuples containing (diameter, count) pairs of descriptions for each of the holes sizes.

Returns:
Vpython:float

Volume [m^3]

Examples

>>> V_multiple_hole_cylinder(0.01, 0.1, [(0.005, 1)])
5.890486225480862e-06