Пример #1
0
def test_bend_rounded_Miller():
    # Miller examples - 9.12
    D = .6
    Re = Reynolds(V=4, D=D, nu=1.14E-6)
    kwargs = dict(Di=D, bend_diameters=2, angle=90,  Re=Re, roughness=.02E-3)

    K = bend_rounded_Miller(L_unimpeded=30*D, **kwargs)
    assert_close(K, 0.1513266131915296, rtol=1e-4)# 0.150 in Miller- 1% difference due to fd
    K = bend_rounded_Miller(L_unimpeded=0*D, **kwargs)
    assert_close(K, 0.1414607344374372, rtol=1e-4) # 0.135 in Miller - Difference mainly from Co interpolation method, OK with that
    K = bend_rounded_Miller(L_unimpeded=2*D, **kwargs)
    assert_close(K, 0.09343184457353562, rtol=1e-4) # 0.093 in miller
Пример #2
0
def calculate(doc_original):
    doc = deepcopy(doc_original)
    treeUnitConvert(doc, doc['units'], SI_UNITS)

    K_pipe = 0
    K_fittings = 0
    K_entry = 0
    K_exit = 0
    K_fixed = 0
    K_LbyD = 0
    K_total = 0
    LbyD_total = 0
    deltaP_fixed = 0
    deltaP_total = 0

    size_definition = doc['input']['pipe']['size_definition']['_val']
    if (size_definition == "NPS"):
        nps = float(doc['input']['pipe']['NPS']['_val'])
        schedule = doc['input']['pipe']['schedule']['_val']
        NPS, Di, Do, t = nearest_pipe(NPS=nps, schedule=schedule)
    else:
        Di = float(doc['input']['pipe']['Dia_inner']['_val'])
        NPS = math.nan
        Do = math.nan
        t = math.nan

    doc['result'].update({'Di': {'_val': str(Di), '_dim': 'length_mili'}})
    doc['result'].update({'Do': {'_val': str(Do), '_dim': 'length_mili'}})
    doc['result'].update({'t': {'_val': str(t), '_dim': 'length_mili'}})

    area = math.pi * pow(Di / 2, 2)
    Q = float(doc['input']['fluidData']['Q']['_val'])
    V = roundit(Q / area)
    doc['result'].update({'V': {'_val': str(V), '_dim': 'speed'}})

    mu = float(doc['input']['fluidData']['mu']['_val'])
    rho = float(doc['input']['fluidData']['rho']['_val'])
    Re = roundit(Reynolds(V=V, D=Di, rho=rho, mu=mu))
    doc['result'].update({'Re': {'_val': str(Re)}})

    Hdyn = roundit(rho * pow(V, 2) / 2)
    doc['result'].update({'Hdyn': {'_val': str(Hdyn), '_dim': 'length'}})

    #K calculation for straigth pipe
    roughness_basis = doc['input']['pipe']['roughness_basis']['_val']
    if (roughness_basis == "Material"):
        material = doc['input']['pipe']['material']['_val']
        roughness = get_roughness(material)
    else:
        roughness = float(doc['input']['pipe']['roughness']['_val'])

    eD = roughness / Di

    doc['result'].update({'eD': {'_val': str(eD)}})
    fd = roundit(friction_factor(Re=Re, eD=eD, Method="Moody"))
    doc['result'].update({'fd_Moody': {'_val': str(fd)}})

    length = float(doc['input']['pipe']['length']['_val'])
    K_pipe = roundit(K_from_f(fd=fd, L=length, D=Di))
    doc['result'].update({'K_pipe': {'_val': str(K_pipe)}})
    deltaP_pipe = roundit(dP_from_K(K_pipe, rho, V))
    doc['result'].update(
        {'deltaP_pipe': {
            '_val': str(deltaP_pipe),
            '_dim': 'pressure'
        }})

    #calculating pressure drop for entrance

    entry_type = doc['input']['entrance']['entry_type']['_val']
    print('entry type is')
    print(entry_type)
    if entry_type == 'none':
        K_entry = 0
    elif entry_type == 'Sharp':
        K_entry = fluids.fittings.entrance_sharp()
    elif entry_type == 'Rounded':
        Rc = float(doc['input']['entrance']['Rc']['_val'])
        K_entry = fluids.fittings.entrance_rounded(Di, Rc)
    elif entry_type == 'Angled':
        angle_radians = float(doc['input']['entrance']['angle']['_val'])
        angle = angle_radians * 57.2958
        K_entry = fluids.fittings.entrance_angled(angle)
    elif entry_type == 'Projecting':
        wall_thickness = float(
            doc['input']['entrance']['wall_thickness']['_val'])
        K_entry = fluids.fittings.entrance_distance(Di, wall_thickness)

    K_entry = roundit(K_entry)
    doc['result'].update({'K_entry': {'_val': str(K_entry)}})
    deltaP_entry = roundit(dP_from_K(K_entry, rho, V))
    doc['result'].update(
        {'deltaP_entry': {
            '_val': str(deltaP_entry),
            '_dim': 'pressure'
        }})

    #calculating pressure drop for exit
    exit_type = doc['input']['exit']['exit_type']['_val']
    print('exit_type')
    print(exit_type)
    if (exit_type == 'Normal'):
        K_exit = exit_normal()
    else:
        K_exit = 0

    K_exit = roundit(K_exit)
    doc['result'].update({'K_exit': {'_val': str(K_exit)}})
    deltaP_exit = roundit(dP_from_K(K_exit, rho, V))
    doc['result'].update(
        {'deltaP_exit': {
            '_val': str(deltaP_exit),
            '_dim': 'pressure'
        }})

    #calculating pressure drop for fittings
    fittings_list = doc['input']['fittings']
    for fitting in fittings_list:
        name = get_hooper_list()[fitting['index']]
        Di_inch = Di * 39.3701
        K_fitting = Hooper2K(Di_inch, Re, name=name)
        K_fittings += K_fitting * fitting['quantity']

    K_fittings = roundit(K_fittings)
    doc['result'].update({'K_fittings': {'_val': str(K_fittings)}})
    deltaP_fittings = dP_from_K(K_fittings, rho, V)
    deltaP_fittings = roundit(deltaP_fittings)
    doc['result'].update({
        'deltaP_fittings': {
            '_val': str(deltaP_fittings),
            '_dim': 'pressure'
        }
    })

    #calculating pressure drop for sharp contractions
    deltaP_contractions_sharp = 0
    contractions_sharp = doc['input']['contractions_sharp']['_list']
    for contraction in contractions_sharp:
        D1 = contraction['D1']
        D2 = contraction['D2']
        A2 = 3.1416 * (D2**2) / 4
        V2 = Q / A2
        K_contraction = fluids.fittings.contraction_sharp(D1, D2)
        deltaP = dP_from_K(K_contraction, rho, V2)
        deltaP_contractions_sharp += deltaP

    deltaP_contractions_sharp = roundit(deltaP_contractions_sharp)
    doc['result'].update({
        'deltaP_contractions_sharp': {
            '_val': str(deltaP_contractions_sharp),
            '_dim': 'pressure'
        }
    })

    #calculating pressure drop for rounded contractions
    deltaP_contractions_rounded = 0
    contractions_rounded = doc['input']['contractions_rounded']['_list']
    for contraction in contractions_rounded:
        D1 = contraction['D1']
        D2 = contraction['D2']
        Rc = contraction['Rc']
        A2 = 3.1416 * (D2**2) / 4
        V2 = Q / A2
        K_contraction = fluids.fittings.contraction_round(D1, D2, Rc)
        deltaP = dP_from_K(K_contraction, rho, V2)
        deltaP_contractions_rounded += deltaP

    deltaP_contractions_rounded = roundit(deltaP_contractions_rounded)
    doc['result'].update({
        'deltaP_contractions_rounded': {
            '_val': str(deltaP_contractions_rounded),
            '_dim': 'pressure'
        }
    })

    #calculating pressure drop for conical contractions
    deltaP_contractions_conical = 0
    contractions_conical = doc['input']['contractions_conical']['_list']
    for contraction in contractions_conical:
        D1 = contraction['D1']
        D2 = contraction['D2']
        L = contraction['L']
        A2 = 3.1416 * (D2**2) / 4
        V2 = Q / A2
        K_contraction = fluids.fittings.contraction_conical(D1, D2, fd=fd, l=L)
        deltaP = dP_from_K(K_contraction, rho, V2)
        deltaP_contractions_conical += deltaP

    deltaP_contractions_conical = roundit(deltaP_contractions_conical)
    doc['result'].update({
        'deltaP_contractions_conical': {
            '_val': str(deltaP_contractions_conical),
            '_dim': 'pressure'
        }
    })

    #calculating pressure drop for pipe reducers contractions
    deltaP_contractions_reducer = 0
    contractions_reducer = doc['input']['contractions_reducer']['_list']
    for contraction in contractions_reducer:
        reducer_size = contraction['reducer_size']
        D1, D2, L = reducer_dimensions(reducer_size)
        A2 = 3.1416 * (D2**2) / 4
        V2 = Q / A2
        K_contraction = fluids.fittings.contraction_conical(D1, D2, fd=fd, l=L)
        deltaP = dP_from_K(K_contraction, rho, V2)
        deltaP_contractions_reducer += deltaP

    deltaP_contractions_reducer = roundit(deltaP_contractions_reducer)
    doc['result'].update({
        'deltaP_contractions_reducer': {
            '_val': str(deltaP_contractions_reducer),
            '_dim': 'pressure'
        }
    })

    # calculating total pressure drop in all contractions
    deltaP_contractions = deltaP_contractions_sharp + deltaP_contractions_rounded + deltaP_contractions_conical + deltaP_contractions_reducer
    deltaP_contractions = roundit(deltaP_contractions)
    doc['result'].update({
        'deltaP_contractions': {
            '_val': str(deltaP_contractions),
            '_dim': 'pressure'
        }
    })

    #calculating pressure drop for sharp expansions
    deltaP_expansions_sharp = 0
    expansions_sharp = doc['input']['expansions_sharp']['_list']
    for contraction in expansions_sharp:
        D1 = contraction['D1']
        D2 = contraction['D2']
        A1 = 3.1416 * (D1**2) / 4
        V1 = Q / A1
        K_contraction = fluids.fittings.diffuser_sharp(D1, D2)
        deltaP = dP_from_K(K_contraction, rho, V1)
        deltaP_expansions_sharp += deltaP

    deltaP_expansions_sharp = roundit(deltaP_expansions_sharp)
    doc['result'].update({
        'deltaP_expansions_sharp': {
            '_val': str(deltaP_expansions_sharp),
            '_dim': 'pressure'
        }
    })

    #calculating pressure drop for conical expansions
    deltaP_expansions_conical = 0
    expansions_conical = doc['input']['expansions_conical']['_list']
    for contraction in expansions_conical:
        D1 = contraction['D1']
        D2 = contraction['D2']
        L = contraction['L']
        A1 = 3.1416 * (D1**2) / 4
        V1 = Q / A1
        K_contraction = fluids.fittings.diffuser_conical(D1, D2, fd=fd, l=L)
        deltaP = dP_from_K(K_contraction, rho, V1)
        deltaP_expansions_conical += deltaP

    deltaP_expansions_conical = roundit(deltaP_expansions_conical)
    doc['result'].update({
        'deltaP_expansions_conical': {
            '_val': str(deltaP_expansions_conical),
            '_dim': 'pressure'
        }
    })

    #calculating pressure drop for pipe reducer expansions
    deltaP_expansions_reducer = 0
    expansions_reducer = doc['input']['expansions_reducer']['_list']
    for contraction in expansions_reducer:
        reducer_size = contraction['reducer_size']
        D2, D1, L = reducer_dimensions(reducer_size)
        A1 = 3.1416 * (D1**2) / 4
        V1 = Q / A1
        K_contraction = fluids.fittings.diffuser_conical(D1, D2, fd=fd, l=L)
        deltaP = dP_from_K(K_contraction, rho, V1)
        deltaP_expansions_reducer += deltaP

    deltaP_expansions_reducer = roundit(deltaP_expansions_reducer)
    doc['result'].update({
        'deltaP_expansions_reducer': {
            '_val': str(deltaP_expansions_reducer),
            '_dim': 'pressure'
        }
    })

    # calculating total pressure drop in all expansions
    deltaP_expansions = deltaP_expansions_sharp + deltaP_expansions_conical + deltaP_expansions_reducer
    doc['result'].update(
        {'deltaP_expansions': {
            '_val': deltaP_expansions,
            '_dim': 'pressure'
        }})

    fixed_K_loss = doc['input']['fixed_K_losses']['_list']
    for loss in fixed_K_loss:
        K_fixed += loss['K'] * loss['quantity']
    deltaP_fixed_K = dP_from_K(K_fixed, rho, V)

    fixed_LbyD_loss = doc['input']['fixed_LbyD_losses']['_list']
    for loss in fixed_LbyD_loss:
        L_D = loss['LbyD']
        K_LbyD += K_from_L_equiv(L_D=L_D, fd=fd) * loss['quantity']
    deltaP_fixed_LbyD = dP_from_K(K_LbyD, rho, V)

    fixed_deltaP_loss = doc['input']['fixed_deltaP_losses']['_list']
    for loss in fixed_deltaP_loss:
        deltaP_fixed += loss['deltaP'] * loss['quantity']
    deltaP_fixed_deltaP = deltaP_fixed

    deltaP_fixed_all = deltaP_fixed_K + deltaP_fixed_LbyD + deltaP_fixed_deltaP

    deltaP_fixed_all = roundit(deltaP_fixed_all)
    doc['result'].update({
        'deltaP_fixed_all': {
            '_val': str(deltaP_fixed_all),
            '_dim': 'pressure'
        }
    })

    deltaP_total = deltaP_pipe + deltaP_entry + deltaP_exit + deltaP_fittings + deltaP_contractions + deltaP_expansions + deltaP_fixed_all
    deltaP_total = roundit(deltaP_total)
    doc['result'].update(
        {'deltaP_total': {
            '_val': str(deltaP_total),
            '_dim': 'pressure'
        }})

    #    doc_original['input'].update(doc['input'])
    doc_original['result'].update(doc['result'])
    treeUnitConvert(doc, SI_UNITS, doc['units'], autoRoundOff=True)
    return True
Пример #3
0
def integrate_drag_sphere(D, rhop, rho, mu, t, V=0, Method=None,
                          distance=False):
    r'''Integrates the velocity and distance traveled by a particle moving
    at a speed which will converge to its terminal velocity.

    Performs an integration of the following expression for acceleration:

    .. math::
        a = \frac{g(\rho_p-\rho_f)}{\rho_p} - \frac{3C_D \rho_f u^2}{4D \rho_p}

    Parameters
    ----------
    D : float
        Diameter of the sphere, [m]
    rhop : float
        Particle density, [kg/m^3]
    rho : float
        Density of the surrounding fluid, [kg/m^3]
    mu : float
        Viscosity of the surrounding fluid [Pa*s]
    t : float
        Time to integrate the particle to, [s]
    V : float
        Initial velocity of the particle, [m/s]
    Method : string, optional
        A string of the function name to use, as in the dictionary
        drag_sphere_correlations
    distance : bool, optional
        Whether or not to calculate the distance traveled and return it as
        well

    Returns
    -------
    v : float
        Velocity of falling sphere after time `t` [m/s]
    x : float, returned only if `distance` == True
        Distance traveled by the falling sphere in time `t`, [m]

    Notes
    -----
    This can be relatively slow as drag correlations can be complex.

    There are analytical solutions available for the Stokes law regime (Re <
    0.3). They were obtained from Wolfram Alpha. [1]_ was not used in the
    derivation, but also describes the derivation fully.

    .. math::
        V(t) = \frac{\exp(-at) (V_0 a + b(\exp(at) - 1))}{a}

    .. math::
        x(t) = \frac{\exp(-a t)\left[V_0 a(\exp(a t) - 1) + b\exp(a t)(a t-1)
        + b\right]}{a^2}

    .. math::
        a = \frac{18\mu_f}{D^2\rho_p}

    .. math::
        b = \frac{g(\rho_p-\rho_f)}{\rho_p}

    The analytical solution will automatically be used if the initial and
    terminal velocity is show the particle's behavior to be laminar. Note
    that this behavior requires that the terminal velocity of the particle be
    solved for - this adds slight (1%) overhead for the cases where particles
    are not laminar.

    Examples
    --------
    >>> integrate_drag_sphere(D=0.001, rhop=2200., rho=1.2, mu=1.78E-5, t=0.5,
    ... V=30, distance=True)
    (9.686465044053, 7.8294546436299)

    References
    ----------
    .. [1] Timmerman, Peter, and Jacobus P. van der Weele. "On the Rise and
       Fall of a Ball with Linear or Quadratic Drag." American Journal of
       Physics 67, no. 6 (June 1999): 538-46. https://doi.org/10.1119/1.19320.
    '''
    # Delayed import of necessaray functions
    from scipy.integrate import odeint, cumtrapz
    import numpy as np
    laminar_initial = Reynolds(V=V, rho=rho, D=D, mu=mu) < 0.01
    v_laminar_end_assumed = v_terminal(D=D, rhop=rhop, rho=rho, mu=mu, Method=Method)
    laminar_end = Reynolds(V=v_laminar_end_assumed, rho=rho, D=D, mu=mu) < 0.01
    if Method == 'Stokes' or (laminar_initial and laminar_end and Method is None):
        try:
            t1 = 18.0*mu/(D*D*rhop)
            t2 = g*(rhop-rho)/rhop
            V_end = exp(-t1*t)*(t1*V + t2*(exp(t1*t) - 1.0))/t1
            x_end = exp(-t1*t)*(V*t1*(exp(t1*t) - 1.0) + t2*exp(t1*t)*(t1*t - 1.0) + t2)/(t1*t1)
            if distance:
                return V_end, x_end
            else:
                return V_end
        except OverflowError:
            # It is only necessary to integrate to terminal velocity
            t_to_terminal = time_v_terminal_Stokes(D, rhop, rho, mu, V0=V, tol=1e-9)
            if t_to_terminal > t:
                raise ValueError('Should never happen')
            V_end, x_end = integrate_drag_sphere(D=D, rhop=rhop, rho=rho, mu=mu, t=t_to_terminal, V=V, Method='Stokes', distance=True)
            # terminal velocity has been reached - V does not change, but x does
            # No reason to believe this isn't working even though it isn't
            # matching the ode solver
            if distance:
                return V_end, x_end + V_end*(t - t_to_terminal)
            else:
                return V_end

            # This is a serious problem for small diameters
            # It would be possible to step slowly, using smaller increments
            # of time to avlid overflows. However, this unfortunately quickly
            # gets much, exponentially, slower than just using odeint because
            # for example solving 10000 seconds might require steps of .0001
            # seconds at a diameter of 1e-7 meters.
#            x = 0.0
#            subdivisions = 10
#            dt = t/subdivisions
#            for i in range(subdivisions):
#                V, dx = integrate_drag_sphere(D=D, rhop=rhop, rho=rho, mu=mu,
#                                              t=dt, V=V, distance=True,
#                                              Method=Method)
#                x += dx
#            if distance:
#                return V, x
#            else:
#                return V

    Re_ish = rho*D/mu
    c1 = g*(rhop-rho)/rhop
    c2 = -0.75*rho/(D*rhop)

    def dv_dt(V, t):
        if V == 0:
            # 64/Re goes to infinity, but gets multiplied by 0 squared.
            t2 = 0.0
        else:
#            t2 = c2*V*V*Stokes(Re_ish*V)
            t2 = c2*V*V*drag_sphere(float(Re_ish*V), Method=Method)
        return c1 + t2

    # Number of intervals for the solution to be solved for; the integrator
    # doesn't care what we give it, but a large number of intervals are needed
    # For an accurate integration of the particle's distance traveled
    pts = 1000 if distance else 2
    ts = np.linspace(0, t, pts)


    # Perform the integration
    Vs = odeint(dv_dt, [V], ts)
    #
    V_end = float(Vs[-1])
    if distance:
        # Calculate the distance traveled
        x = float(cumtrapz(np.ravel(Vs), ts)[-1])
        return V_end, x
    else:
        return V_end
Пример #4
0
def v_terminal(D, rhop, rho, mu, Method=None):
    r'''Calculates terminal velocity of a falling sphere using any drag
    coefficient method supported by `drag_sphere`. The laminar solution for
    Re < 0.01 is first tried; if the resulting terminal velocity does not
    put it in the laminar regime, a numerical solution is used.

    .. math::
        v_t = \sqrt{\frac{4 g d_p (\rho_p-\rho_f)}{3 C_D \rho_f }}

    Parameters
    ----------
    D : float
        Diameter of the sphere, [m]
    rhop : float
        Particle density, [kg/m^3]
    rho : float
        Density of the surrounding fluid, [kg/m^3]
    mu : float
        Viscosity of the surrounding fluid [Pa*s]
    Method : string, optional
        A string of the function name to use, as in the dictionary
        drag_sphere_correlations

    Returns
    -------
    v_t : float
        Terminal velocity of falling sphere [m/s]

    Notes
    -----
    As there are no correlations implemented for Re > 1E6, an error will be
    raised if the numerical solver seeks a solution above that limit.

    The laminar solution is given in [1]_ and is:

    .. math::
        v_t = \frac{g d_p^2 (\rho_p - \rho_f)}{18 \mu_f}

    Examples
    --------
    >>> v_terminal(D=70E-6, rhop=2600., rho=1000., mu=1E-3)
    0.004142497244531304

    Example 7-1 in GPSA handbook, 13th edition:

    >>> from scipy.constants import *
    >>> v_terminal(D=150E-6, rhop=31.2*lb/foot**3, rho=2.07*lb/foot**3,  mu=1.2e-05)/foot
    0.4491992020345101

    The answer reported there is 0.46 ft/sec.

    References
    ----------
    .. [1] Green, Don, and Robert Perry. Perry's Chemical Engineers' Handbook,
       Eighth Edition. McGraw-Hill Professional, 2007.
    .. [2] Rushton, Albert, Anthony S. Ward, and Richard G. Holdich.
       Solid-Liquid Filtration and Separation Technology. 1st edition. Weinheim ;
       New York: Wiley-VCH, 1996.
    '''
    '''The following would be the ideal implementation. The actual function is
    optimized for speed, not readability
    def err(V):
        Re = rho*V*D/mu
        Cd = Barati_high(Re)
        V2 = (4/3.*g*D*(rhop-rho)/rho/Cd)**0.5
        return (V-V2)
    return fsolve(err, 1.)'''
    v_lam = g*D*D*(rhop-rho)/(18*mu)
    Re_lam = Reynolds(V=v_lam, D=D, rho=rho, mu=mu)
    if Re_lam < 0.01 or Method == 'Stokes':
        return v_lam

    Re_almost = rho*D/mu
    main = 4/3.*g*D*(rhop-rho)/rho
    V_max = 1E6/rho/D*mu  # where the correlation breaks down, Re=1E6

    # Begin the solver with 1/100 th the velocity possible at the maximum
    # Reynolds number the correlation is good for
    return secant(_v_terminal_err, V_max*1e-2, xtol=1E-12, args=(Method, Re_almost, main))
Пример #5
0
def h_Ganguli_VDI(m, A, A_min, A_increase, A_fin, A_tube_showing,
                  tube_diameter, fin_diameter, fin_thickness, bare_length,
                  pitch_parallel, pitch_normal, tube_rows, rho, Cp, mu, k,
                  k_fin):
    r'''Calculates the air side heat transfer coefficient for an air cooler
    or other finned tube bundle with the formulas of [1]_ as modified in [2]_.
    
    Inline:
        
    .. math::
        Nu_d = 0.22Re_d^{0.6}\left(\frac{A}{A_{tube,only}}\right)^{-0.15}Pr^{1/3}
        
    Staggered:
        
    .. math::
        Nu_d = 0.38 Re_d^{0.6}\left(\frac{A}{A_{tube,only}}\right)^{-0.15}Pr^{1/3}
            
    Parameters
    ----------
    m : float
        Mass flow rate of air across the tube bank, [kg/s]
    A : float
        Surface area of combined finned and non-finned area exposed for heat
        transfer, [m^2]
    A_min : float
        Minimum air flow area, [m^2]
    A_increase : float
        Ratio of actual surface area to bare tube surface area
        :math:`A_{increase} = \frac{A_{tube}}{A_{bare, total/tube}}`, [-]
    A_fin : float
        Surface area of all fins in the bundle, [m^2]
    A_tube_showing : float
        Area of the bare tube which is exposed in the bundle, [m^2]
    tube_diameter : float
        Diameter of the bare tube, [m]
    fin_diameter : float
        Outer diameter of each tube after including the fin on both sides,
        [m]
    fin_thickness : float
        Thickness of the fins, [m]
    bare_length : float
        Length of bare tube between two fins 
        :math:`\text{bare length} = \text{fin interval} - t_{fin}`, [m]
    pitch_parallel : float
        Distance between tube center along a line parallel to the flow;
        has been called `longitudinal` pitch, `pp`, `s2`, `SL`, and `p2`, [m]
    pitch_normal : float
        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]
    tube_rows : int
        Number of tube rows per bundle, [-]
    rho : float
        Average (bulk) density of air across the tube bank, [kg/m^3]
    Cp : float
        Average (bulk) heat capacity of air across the tube bank, [J/kg/K]
    mu : float
        Average (bulk) viscosity of air across the tube bank, [Pa*s]
    k : float
        Average (bulk) thermal conductivity of air across the tube bank, 
        [W/m/K]
    k_fin : float
        Thermal conductivity of the fin, [W/m/K]
        
    Returns
    -------
    h_bare_tube_basis : float
        Air side heat transfer coefficient on a bare-tube surface area as if 
        there were no fins present basis, [W/K/m^2]

    Notes
    -----
    The VDI modifications were developed in comparison with HTFS and HTRI data
    according to [2]_.
    
    For cases where the tube row count is less than four, the coefficients are
    modified in [2]_. For the inline case, 0.2 replaces 0.22. For the stagered
    cases, the coefficient is 0.2, 0.33, 0.36 for 1, 2, or 3 tube rows 
    respectively.
    
    The model is also showin in [4]_.
    
    Examples
    --------
    Example 12.1 in [3]_:
    
    >>> AC = AirCooledExchanger(tube_rows=4, tube_passes=4, tubes_per_row=56, tube_length=36*foot, 
    ... tube_diameter=1*inch, fin_thickness=0.013*inch, fin_density=10/inch,
    ... angle=30, pitch_normal=2.5*inch, fin_height=0.625*inch, corbels=True)

    >>> h_Ganguli_VDI(m=130.70315, A=AC.A, A_min=AC.A_min, A_increase=AC.A_increase, A_fin=AC.A_fin,
    ... A_tube_showing=AC.A_tube_showing, tube_diameter=AC.tube_diameter,
    ... fin_diameter=AC.fin_diameter, bare_length=AC.bare_length,
    ... fin_thickness=AC.fin_thickness, tube_rows=AC.tube_rows,
    ... pitch_parallel=AC.pitch_parallel, pitch_normal=AC.pitch_normal,
    ... rho=1.2013848, Cp=1009.0188, mu=1.9304793e-05, k=0.027864828, k_fin=238)
    969.2850818578595
    
    References
    ----------
    .. [1] Ganguli, A., S. S. Tung, and J. Taborek. "Parametric Study of
       Air-Cooled Heat Exchanger Finned Tube Geometry." In AIChE Symposium 
       Series, 81:122-28, 1985.
    .. [2] Gesellschaft, V. D. I., ed. VDI Heat Atlas. 2nd edition.
       Berlin; New York:: Springer, 2010.
    .. [3] Serth, Robert W., and Thomas Lestina. Process Heat Transfer: 
       Principles, Applications and Rules of Thumb. Academic Press, 2014.
    .. [4] Kroger, Detlev. Air-Cooled Heat Exchangers and Cooling Towers: 
       Thermal-Flow Performance Evaluation and Design, Vol. 1. Tulsa, Okl:
       PennWell Corp., 2004.
    '''
    V_max = m / (A_min * rho)

    Re = Reynolds(V=V_max, D=tube_diameter, rho=rho, mu=mu)
    Pr = Prandtl(Cp=Cp, mu=mu, k=k)

    if abs(1 - pitch_normal / pitch_parallel
           ) < 0.05:  # in-line, with a tolerance of 0.05 proximity
        if tube_rows < 4:
            coeff = 0.2
        else:
            coeff = 0.22
    else:  # staggered
        if tube_rows == 1:
            coeff = 0.2
        elif tube_rows == 2:
            coeff = 0.33
        elif tube_rows == 3:
            coeff = 0.36
        else:
            coeff = 0.38

    # VDI example shows the ratio is of the total area, to the original bare tube area
    # Serth example would match Nu = 47.22 except for lazy rounding
    Nu = coeff * Re**0.6 * Pr**(1 / 3.) * (A_increase)**-0.15
    h = k / tube_diameter * Nu
    efficiency = fin_efficiency_Kern_Kraus(Do=tube_diameter,
                                           D_fin=fin_diameter,
                                           t_fin=fin_thickness,
                                           k_fin=k_fin,
                                           h=h)
    h_total_area_basis = (efficiency * A_fin + A_tube_showing) / A * h
    h_bare_tube_basis = h_total_area_basis * A_increase
    return h_bare_tube_basis
Пример #6
0
def h_ESDU_low_fin(m,
                   A,
                   A_min,
                   A_increase,
                   A_fin,
                   A_tube_showing,
                   tube_diameter,
                   fin_diameter,
                   fin_thickness,
                   bare_length,
                   pitch_parallel,
                   pitch_normal,
                   tube_rows,
                   rho,
                   Cp,
                   mu,
                   k,
                   k_fin,
                   Pr_wall=None):
    r'''Calculates the air side heat transfer coefficient for an air cooler
    or other finned tube bundle with low fins using the formulas of [1]_ as
    presented in [2]_ (and also [3]_).
    
    .. math::
        Nu = 0.183Re^{0.7} \left(\frac{\text{bare length}}{\text{fin height}}
        \right)^{0.36}\left(\frac{p_1}{D_{o}}\right)^{0.06}
        \left(\frac{\text{fin height}}{D_o}\right)^{0.11}
        Pr^{0.36} \cdot F_1\cdot F_2
        
    .. math::
        h_{A,total} = \frac{\eta A_{fin} + A_{bare, showing}}{A_{total}} h
        
    .. math::
        h_{bare,total} = A_{increase} h_{A,total}
        
    Parameters
    ----------
    m : float
        Mass flow rate of air across the tube bank, [kg/s]
    A : float
        Surface area of combined finned and non-finned area exposed for heat
        transfer, [m^2]
    A_min : float
        Minimum air flow area, [m^2]
    A_increase : float
        Ratio of actual surface area to bare tube surface area
        :math:`A_{increase} = \frac{A_{tube}}{A_{bare, total/tube}}`, [-]
    A_fin : float
        Surface area of all fins in the bundle, [m^2]
    A_tube_showing : float
        Area of the bare tube which is exposed in the bundle, [m^2]
    tube_diameter : float
        Diameter of the bare tube, [m]
    fin_diameter : float
        Outer diameter of each tube after including the fin on both sides,
        [m]
    fin_thickness : float
        Thickness of the fins, [m]
    bare_length : float
        Length of bare tube between two fins 
        :math:`\text{bare length} = \text{fin interval} - t_{fin}`, [m]
    pitch_parallel : float
        Distance between tube center along a line parallel to the flow;
        has been called `longitudinal` pitch, `pp`, `s2`, `SL`, and `p2`, [m]
    pitch_normal : float
        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]
    tube_rows : int
        Number of tube rows per bundle, [-]
    rho : float
        Average (bulk) density of air across the tube bank, [kg/m^3]
    Cp : float
        Average (bulk) heat capacity of air across the tube bank, [J/kg/K]
    mu : float
        Average (bulk) viscosity of air across the tube bank, [Pa*s]
    k : float
        Average (bulk) thermal conductivity of air across the tube bank, 
        [W/m/K]
    k_fin : float
        Thermal conductivity of the fin, [W/m/K]
    Pr_wall : float, optional
        Prandtl number at the wall temperature; provide if a correction with  
        the defaults parameters is desired; otherwise apply the correction
        elsewhere, [-]
        
    Returns
    -------
    h_bare_tube_basis : float
        Air side heat transfer coefficient on a bare-tube surface area as if 
        there were no fins present basis, [W/K/m^2]

    Notes
    -----
    The tube-row count correction factor `F2` can be disabled by setting `tube_rows`
    to 10. The property correction factor `F1` can be disabled by not specifying
    `Pr_wall`. A Prandtl number exponent of 0.26 is recommended in [1]_ for 
    heating and cooling for both liquids and gases.
    
    There is a third correction factor in [1]_ for tube angles not 30, 45, or
    60 degrees, but it is not fully explained and it is not shown in [2]_. 
    Another correction factor is in [2]_ for flow at an angle; however it would
    not make sense to apply it to finned tube banks due to the blockage by the
    fins.
    
    Examples
    --------
    >>> AC = AirCooledExchanger(tube_rows=4, tube_passes=4, tubes_per_row=8, tube_length=0.5, 
    ... tube_diameter=0.0164, fin_thickness=0.001, fin_density=1/0.003,
    ... pitch_normal=0.0313, pitch_parallel=0.0271, fin_height=0.0041, corbels=True)
    
    >>> h_ESDU_low_fin(m=0.914, A=AC.A, A_min=AC.A_min, A_increase=AC.A_increase, A_fin=AC.A_fin,
    ... A_tube_showing=AC.A_tube_showing, tube_diameter=AC.tube_diameter,
    ... fin_diameter=AC.fin_diameter, bare_length=AC.bare_length,
    ... fin_thickness=AC.fin_thickness, tube_rows=AC.tube_rows,
    ... pitch_normal=AC.pitch_normal, pitch_parallel=AC.pitch_parallel, 
    ... rho=1.217, Cp=1007., mu=1.8E-5, k=0.0253, k_fin=15)
    553.853836470948

    References
    ----------
    .. [1] Hewitt, G. L. Shires, T. Reg Bott G. F., George L. Shires, and T.
       R. Bott. Process Heat Transfer. 1st edition. Boca Raton: CRC Press, 
       1994.
    .. [2] "High-Fin Staggered Tube Banks: Heat Transfer and Pressure Drop for
       Turbulent Single Phase Gas Flow." ESDU 86022 (October 1, 1986). 
    .. [3] Rabas, T. J., and J. Taborek. "Survey of Turbulent Forced-Convection
       Heat Transfer and Pressure Drop Characteristics of Low-Finned Tube Banks
       in Cross Flow."  Heat Transfer Engineering 8, no. 2 (January 1987): 
       49-62.
    '''
    fin_height = 0.5 * (fin_diameter - tube_diameter)

    V_max = m / (A_min * rho)
    Re = Reynolds(V=V_max, D=tube_diameter, rho=rho, mu=mu)
    Pr = Prandtl(Cp=Cp, mu=mu, k=k)
    Nu = (0.183 * Re**0.7 * (bare_length / fin_height)**0.36 *
          (pitch_normal / fin_diameter)**0.06 *
          (fin_height / fin_diameter)**0.11 * Pr**0.36)

    staggered = abs(1 - pitch_normal / pitch_parallel) > 0.05
    F2 = ESDU_tube_row_correction(tube_rows=tube_rows, staggered=staggered)
    Nu *= F2
    if Pr_wall is not None:
        F1 = wall_factor(Pr=Pr,
                         Pr_wall=Pr_wall,
                         Pr_heating_coeff=0.26,
                         Pr_cooling_coeff=0.26,
                         property_option=WALL_FACTOR_PRANDTL)
        Nu *= F1

    h = k / tube_diameter * Nu
    efficiency = fin_efficiency_Kern_Kraus(Do=tube_diameter,
                                           D_fin=fin_diameter,
                                           t_fin=fin_thickness,
                                           k_fin=k_fin,
                                           h=h)
    h_total_area_basis = (efficiency * A_fin + A_tube_showing) / A * h
    h_bare_tube_basis = h_total_area_basis * A_increase
    return h_bare_tube_basis
Пример #7
0
def h_ESDU_high_fin(m,
                    A,
                    A_min,
                    A_increase,
                    A_fin,
                    A_tube_showing,
                    tube_diameter,
                    fin_diameter,
                    fin_thickness,
                    bare_length,
                    pitch_parallel,
                    pitch_normal,
                    tube_rows,
                    rho,
                    Cp,
                    mu,
                    k,
                    k_fin,
                    Pr_wall=None):
    r'''Calculates the air side heat transfer coefficient for an air cooler
    or other finned tube bundle with the formulas of [2]_ as presented in [1]_.

    .. math::
        Nu = 0.242 Re^{0.658} \left(\frac{\text{bare length}}
        {\text{fin height}}\right)^{0.297}
        \left(\frac{P_1}{P_2}\right)^{-0.091} P_r^{1/3}\cdot F_1\cdot F_2
                
    .. math::
        h_{A,total} = \frac{\eta A_{fin} + A_{bare, showing}}{A_{total}} h
        
    .. math::
        h_{bare,total} = A_{increase} h_{A,total}

    Parameters
    ----------
    m : float
        Mass flow rate of air across the tube bank, [kg/s]
    A : float
        Surface area of combined finned and non-finned area exposed for heat
        transfer, [m^2]
    A_min : float
        Minimum air flow area, [m^2]
    A_increase : float
        Ratio of actual surface area to bare tube surface area
        :math:`A_{increase} = \frac{A_{tube}}{A_{bare, total/tube}}`, [-]
    A_fin : float
        Surface area of all fins in the bundle, [m^2]
    A_tube_showing : float
        Area of the bare tube which is exposed in the bundle, [m^2]
    tube_diameter : float
        Diameter of the bare tube, [m]
    fin_diameter : float
        Outer diameter of each tube after including the fin on both sides,
        [m]
    fin_thickness : float
        Thickness of the fins, [m]
    bare_length : float
        Length of bare tube between two fins 
        :math:`\text{bare length} = \text{fin interval} - t_{fin}`, [m]
    pitch_parallel : float
        Distance between tube center along a line parallel to the flow;
        has been called `longitudinal` pitch, `pp`, `s2`, `SL`, and `p2`, [m]
    pitch_normal : float
        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]
    tube_rows : int
        Number of tube rows per bundle, [-]
    rho : float
        Average (bulk) density of air across the tube bank, [kg/m^3]
    Cp : float
        Average (bulk) heat capacity of air across the tube bank, [J/kg/K]
    mu : float
        Average (bulk) viscosity of air across the tube bank, [Pa*s]
    k : float
        Average (bulk) thermal conductivity of air across the tube bank, 
        [W/m/K]
    k_fin : float
        Thermal conductivity of the fin, [W/m/K]
    Pr_wall : float, optional
        Prandtl number at the wall temperature; provide if a correction with  
        the defaults parameters is desired; otherwise apply the correction
        elsewhere, [-]

    Returns
    -------
    h_bare_tube_basis : float
        Air side heat transfer coefficient on a bare-tube surface area as if 
        there were no fins present basis, [W/K/m^2]

    Notes
    -----
    The tube-row count correction factor is 1 for four or more rows, 0.92 for
    three rows, 0.84 for two rows, and 0.76 for one row according to [1]_.
    
    The property correction factor can be disabled by not specifying
    `Pr_wall`. A Prandtl number exponent of 0.26 is recommended in [1]_ for 
    heating and cooling for both liquids and gases.
    
    Examples
    --------
    >>> AC = AirCooledExchanger(tube_rows=4, tube_passes=4, tubes_per_row=20, tube_length=3, 
    ... tube_diameter=1*inch, fin_thickness=0.000406, fin_density=1/0.002309,
    ... pitch_normal=.06033, pitch_parallel=.05207,
    ... fin_height=0.0159, tube_thickness=(.0254-.0186)/2,
    ... bundles_per_bay=1, parallel_bays=1, corbels=True)
    
    >>> h_ESDU_high_fin(m=21.56, A=AC.A, A_min=AC.A_min, A_increase=AC.A_increase, A_fin=AC.A_fin,
    ... A_tube_showing=AC.A_tube_showing, tube_diameter=AC.tube_diameter,
    ... fin_diameter=AC.fin_diameter, bare_length=AC.bare_length,
    ... fin_thickness=AC.fin_thickness, tube_rows=AC.tube_rows,
    ... pitch_normal=AC.pitch_normal, pitch_parallel=AC.pitch_parallel, 
    ... rho=1.161, Cp=1007., mu=1.85E-5, k=0.0263, k_fin=205)
    1390.888918049757

    References
    ----------
    .. [1] Hewitt, G. L. Shires, T. Reg Bott G. F., George L. Shires, and T.
       R. Bott. Process Heat Transfer. 1st edition. Boca Raton: CRC Press, 
       1994.
    .. [2] "High-Fin Staggered Tube Banks: Heat Transfer and Pressure Drop for
       Turbulent Single Phase Gas Flow." ESDU 86022 (October 1, 1986). 
    .. [3] Rabas, T. J., and J. Taborek. "Survey of Turbulent Forced-Convection
       Heat Transfer and Pressure Drop Characteristics of Low-Finned Tube Banks
       in Cross Flow."  Heat Transfer Engineering 8, no. 2 (January 1987): 
       49-62.
    '''
    fin_height = 0.5 * (fin_diameter - tube_diameter)

    V_max = m / (A_min * rho)
    Re = Reynolds(V=V_max, D=tube_diameter, rho=rho, mu=mu)
    Pr = Prandtl(Cp=Cp, mu=mu, k=k)
    Nu = 0.242 * Re**0.658 * (bare_length / fin_height)**0.297 * (
        pitch_normal / pitch_parallel)**-0.091 * Pr**(1 / 3.)

    if tube_rows < 2:
        F2 = 0.76
    elif tube_rows < 3:
        F2 = 0.84
    elif tube_rows < 4:
        F2 = 0.92
    else:
        F2 = 1.0

    Nu *= F2
    if Pr_wall is not None:
        F1 = wall_factor(Pr=Pr,
                         Pr_wall=Pr_wall,
                         Pr_heating_coeff=0.26,
                         Pr_cooling_coeff=0.26,
                         property_option=WALL_FACTOR_PRANDTL)
        Nu *= F1

    h = k / tube_diameter * Nu
    efficiency = fin_efficiency_Kern_Kraus(Do=tube_diameter,
                                           D_fin=fin_diameter,
                                           t_fin=fin_thickness,
                                           k_fin=k_fin,
                                           h=h)

    h_total_area_basis = (efficiency * A_fin + A_tube_showing) / A * h
    h_bare_tube_basis = h_total_area_basis * A_increase
    return h_bare_tube_basis
Пример #8
0
def h_Briggs_Young(m, A, A_min, A_increase, A_fin, A_tube_showing,
                   tube_diameter, fin_diameter, fin_thickness, bare_length,
                   rho, Cp, mu, k, k_fin):
    r'''Calculates the air side heat transfer coefficient for an air cooler
    or other finned tube bundle with the formulas of Briggs and Young [1], [2]_
    [3]_.
    
    .. math::
        Nu = 0.134Re^{0.681} Pr^{0.33}\left(\frac{S}{h}\right)^{0.2}
        \left(\frac{S}{b}\right)^{0.1134}
        
    Parameters
    ----------
    m : float
        Mass flow rate of air across the tube bank, [kg/s]
    A : float
        Surface area of combined finned and non-finned area exposed for heat
        transfer, [m^2]
    A_min : float
        Minimum air flow area, [m^2]
    A_increase : float
        Ratio of actual surface area to bare tube surface area
        :math:`A_{increase} = \frac{A_{tube}}{A_{bare, total/tube}}`, [-]
    A_fin : float
        Surface area of all fins in the bundle, [m^2]
    A_tube_showing : float
        Area of the bare tube which is exposed in the bundle, [m^2]
    tube_diameter : float
        Diameter of the bare tube, [m]
    fin_diameter : float
        Outer diameter of each tube after including the fin on both sides,
        [m]
    fin_thickness : float
        Thickness of the fins, [m]
    bare_length : float
        Length of bare tube between two fins 
        :math:`\text{bare length} = \text{fin interval} - t_{fin}`, [m]
    rho : float
        Average (bulk) density of air across the tube bank, [kg/m^3]
    Cp : float
        Average (bulk) heat capacity of air across the tube bank, [J/kg/K]
    mu : float
        Average (bulk) viscosity of air across the tube bank, [Pa*s]
    k : float
        Average (bulk) thermal conductivity of air across the tube bank, 
        [W/m/K]
    k_fin : float
        Thermal conductivity of the fin, [W/m/K]
        
    Returns
    -------
    h_bare_tube_basis : float
        Air side heat transfer coefficient on a bare-tube surface area as if 
        there were no fins present basis, [W/K/m^2]

    Notes
    -----
    The limits on this equation are :math:`1000 < Re < `8000`, 
    11.13 mm :math:`< D_o < ` 40.89 mm, 1.42 mm < fin height < 16.57 mm,
    0.33 mm < fin thickness < 2.02 mm, 1.30 mm < fin pitch < 4.06 mm, and 
    24.49 mm < normal pitch < 111 mm.
    
    Examples
    --------
    >>> AC = AirCooledExchanger(tube_rows=4, tube_passes=4, tubes_per_row=20, tube_length=3, 
    ... tube_diameter=1*inch, fin_thickness=0.000406, fin_density=1/0.002309,
    ... pitch_normal=.06033, pitch_parallel=.05207,
    ... fin_height=0.0159, tube_thickness=(.0254-.0186)/2,
    ... bundles_per_bay=1, parallel_bays=1, corbels=True)
    
    >>> h_Briggs_Young(m=21.56, A=AC.A, A_min=AC.A_min, A_increase=AC.A_increase, A_fin=AC.A_fin,
    ... A_tube_showing=AC.A_tube_showing, tube_diameter=AC.tube_diameter,
    ... fin_diameter=AC.fin_diameter, bare_length=AC.bare_length,
    ... fin_thickness=AC.fin_thickness,
    ... rho=1.161, Cp=1007., mu=1.85E-5, k=0.0263, k_fin=205)
    1422.8722403237716

    References
    ----------
    .. [1] Briggs, D.E., and Young, E.H., 1963, "Convection Heat Transfer and
       Pressure Drop of Air Flowing across Triangular Banks of Finned Tubes",
       Chemical Engineering Progress Symp., Series 41, No. 59. Chem. Eng. Prog.
       Symp. Series No. 41, "Heat Transfer - Houston".
    .. [2] Mukherjee, R., and Geoffrey Hewitt. Practical Thermal Design of 
       Air-Cooled Heat Exchangers. New York: Begell House Publishers Inc.,U.S.,
       2007.
    .. [3] Kroger, Detlev. Air-Cooled Heat Exchangers and Cooling Towers: 
       Thermal-Flow Performance Evaluation and Design, Vol. 1. Tulsa, Okl:
       PennWell Corp., 2004.
    '''
    fin_height = 0.5 * (fin_diameter - tube_diameter)

    V_max = m / (A_min * rho)

    Re = Reynolds(V=V_max, D=tube_diameter, rho=rho, mu=mu)
    Pr = Prandtl(Cp=Cp, mu=mu, k=k)

    Nu = 0.134 * Re**0.681 * Pr**(1 / 3.) * (bare_length / fin_height)**0.2 * (
        bare_length / fin_thickness)**0.1134

    h = k / tube_diameter * Nu
    efficiency = fin_efficiency_Kern_Kraus(Do=tube_diameter,
                                           D_fin=fin_diameter,
                                           t_fin=fin_thickness,
                                           k_fin=k_fin,
                                           h=h)
    h_total_area_basis = (efficiency * A_fin + A_tube_showing) / A * h
    h_bare_tube_basis = h_total_area_basis * A_increase

    return h_bare_tube_basis
Пример #9
0
def Aggour(m, x, alpha, D, rhol, Cpl, kl, mu_b, mu_w=None, L=None,
           turbulent=None):
    r'''Calculates the two-phase non-boiling laminar heat transfer coefficient
    of a liquid and gas flowing inside a tube of any inclination, as in [1]_
    and reviewed in [2]_.

    Laminar for Rel <= 2000:

    .. math::
        h_{TP} = 1.615\frac{k_l}{D}\left(\frac{Re_l Pr_l D}{L}\right)^{1/3}
        \left(\frac{\mu_b}{\mu_w}\right)^{0.14}

    Turbulent for Rel > 2000:

    .. math::
        h_{TP} = 0.0155\frac{k_l}{D} Pr_l^{0.5} Re_l^{0.83}

    .. math::
        Re_l = \frac{\rho_l v_l D}{\mu_l}

    .. math::
        V_l = \frac{V_{ls}}{1-\alpha}

    Parameters
    ----------
    m : float
        Mass flow rate [kg/s]
    x : float
        Quality at the specific tube interval [-]
    alpha : float
        Void fraction in the tube, [-]
    D : float
        Diameter of the tube [m]
    rhol : float
        Density of the liquid [kg/m^3]
    Cpl : float
        Constant-pressure heat capacity of liquid [J/kg/K]
    kl : float
        Thermal conductivity of liquid [W/m/K]
    mu_b : float
        Viscosity of liquid at bulk conditions (average of inlet/outlet
        temperature) [Pa*s]
    mu_w : float, optional
        Viscosity of liquid at wall temperature [Pa*s]
    L : float, optional
        Length of the tube, [m]
    turbulent : bool or None, optional
        Whether or not to force the correlation to return the turbulent
        result; will return the laminar regime if False

    Returns
    -------
    h : float
        Heat transfer coefficient [W/m^2/K]

    Notes
    -----
    Developed with mixtures of air-water, helium-water, and freon-12-water and
    vertical tests. Studied flow patterns were bubbly, slug, annular,
    bubbly-slug, and slug-annular regimes. Superficial velocity ratios ranged
    from 0.02 to 470.

    A viscosity correction is only suggested for the laminar regime.
    If the viscosity at the wall temperature is not given, the liquid viscosity
    correction is not applied.

    Examples
    --------
    >>> Aggour(m=1, x=.9, D=.3, alpha=.9, rhol=1000, Cpl=2300, kl=.6, mu_b=1E-3)
    420.9347146885667

    References
    ----------
    .. [1] Aggour, Mohamed A. Hydrodynamics and Heat Transfer in Two-Phase
       Two-Component Flows, Ph.D. Thesis, University of Manutoba, Canada
       (1978). http://mspace.lib.umanitoba.ca/xmlui/handle/1993/14171.
    .. [2] Dongwoo Kim, Venkata K. Ryali, Afshin J. Ghajar, Ronald L.
       Dougherty. "Comparison of 20 Two-Phase Heat Transfer Correlations with
       Seven Sets of Experimental Data, Including Flow Pattern and Tube
       Inclination Effects." Heat Transfer Engineering 20, no. 1 (February 1,
       1999): 15-40. doi:10.1080/014576399271691.
    '''
    Vls = m*(1-x)/(rhol*pi/4*D**2)
    Vl = Vls/(1.-alpha)

    Prl = Prandtl(Cp=Cpl, k=kl, mu=mu_b)
    Rel = Reynolds(V=Vl, D=D, rho=rhol, mu=mu_b)

    if turbulent or (Rel > 2000.0 and turbulent is None):
        hl = 0.0155*(kl/D)*Rel**0.83*Prl**0.5
        return hl*(1-alpha)**-0.83
    else:
        hl = 1.615*(kl/D)*(Rel*Prl*D/L)**(1/3.)
        if mu_w:
            hl *= (mu_b/mu_w)**0.14
        return hl*(1.0 - alpha)**(-1/3.)
Пример #10
0
def Martin_Sims(m, x, D, rhol, rhog, hl=None,
                Cpl=None, kl=None, mu_b=None, mu_w=None, L=None):
    r'''Calculates the two-phase non-boiling heat transfer coefficient of a
    liquid and gas flowing inside a tube of any inclination, as in [1]_ and
    reviewed in [2]_.

    .. math::
        \frac{h_{TP}}{h_l} = 1 + 0.64\sqrt{\frac{V_{gs}}{V_{ls}}}

    Parameters
    ----------
    m : float
        Mass flow rate [kg/s]
    x : float
        Quality at the specific tube interval []
    D : float
        Diameter of the tube [m]
    rhol : float
        Density of the liquid [kg/m^3]
    rhog : float
        Density of the gas [kg/m^3]
    hl : float, optional
        Liquid-phase heat transfer coefficient as described below, [W/m^2/K]
    Cpl : float, optional
        Constant-pressure heat capacity of liquid [J/kg/K]
    kl : float, optional
        Thermal conductivity of liquid [W/m/K]
    mu_b : float, optional
        Viscosity of liquid at bulk conditions (average of inlet/outlet
        temperature) [Pa*s]
    mu_w : float, optional
        Viscosity of liquid at wall temperature [Pa*s]
    L : float, optional
        Length of the tube [m]

    Returns
    -------
    h : float
        Heat transfer coefficient [W/m^2/K]

    Notes
    -----
    No specific suggestion for how to calculate the liquid-phase heat transfer
    coefficient is given in [1]_; [2]_ suggests to use the same procedure as
    in `Knott`.

    Examples
    --------
    >>> Martin_Sims(m=1, x=.9, D=.3, rhol=1000, rhog=2.5, hl=141.2)
    5563.280000000001
    >>> Martin_Sims(m=1, x=.9, D=.3, rhol=1000, rhog=2.5, Cpl=2300, kl=.6,
    ... mu_b=1E-3, mu_w=1.2E-3, L=24)
    5977.505465781747

    References
    ----------
    .. [1] Martin, B. W, and G. E Sims. "Forced Convection Heat Transfer to
       Water with Air Injection in a Rectangular Duct." International Journal
       of Heat and Mass Transfer 14, no. 8 (August 1, 1971): 1115-34.
       doi:10.1016/0017-9310(71)90208-0.
    .. [2] Dongwoo Kim, Venkata K. Ryali, Afshin J. Ghajar, Ronald L.
       Dougherty. "Comparison of 20 Two-Phase Heat Transfer Correlations with
       Seven Sets of Experimental Data, Including Flow Pattern and Tube
       Inclination Effects." Heat Transfer Engineering 20, no. 1 (February 1,
       1999): 15-40. doi:10.1080/014576399271691.
    '''
    Vgs = m*x/(rhog*pi/4*D**2)
    Vls = m*(1-x)/(rhol*pi/4*D**2)
    if hl is None:
        V = Vgs + Vls # Net velocity
        Re = Reynolds(V=V, D=D, rho=rhol, mu=mu_b)
        Pr = Prandtl(Cp=Cpl, k=kl, mu=mu_b)
        Nul = laminar_entry_Seider_Tate(Re=Re, Pr=Pr, L=L, Di=D, mu=mu_b, mu_w=mu_w)
        hl = Nul*kl/D
    return hl*(1.0 + 0.64*(Vgs/Vls)**0.5)
Пример #11
0
def Knott(m, x, D, rhol, rhog, Cpl=None, kl=None, mu_b=None, mu_w=None, L=None,
          hl=None):
    r'''Calculates the two-phase non-boiling heat transfer coefficient of a
    liquid and gas flowing inside a tube of any inclination, as in [1]_ and
    reviewed in [2]_.

    Either a specified `hl` is required, or `Cpl`, `kl`, `mu_b`, `mu_w` and
    `L` are required to calculate `hl`.

    .. math::
        \frac{h_{TP}}{h_l} = \left(1 + \frac{V_{gs}}{V_{ls}}\right)^{1/3}

    Parameters
    ----------
    m : float
        Mass flow rate [kg/s]
    x : float
        Quality at the specific tube interval [-]
    D : float
        Diameter of the tube [m]
    rhol : float
        Density of the liquid [kg/m^3]
    rhog : float
        Density of the gas [kg/m^3]
    Cpl : float, optional
        Constant-pressure heat capacity of liquid [J/kg/K]
    kl : float, optional
        Thermal conductivity of liquid [W/m/K]
    mu_b : float, optional
        Viscosity of liquid at bulk conditions (average of inlet/outlet
        temperature) [Pa*s]
    mu_w : float, optional
        Viscosity of liquid at wall temperature [Pa*s]
    L : float, optional
        Length of the tube [m]
    hl : float, optional
        Liquid-phase heat transfer coefficient as described below, [W/m^2/K]

    Returns
    -------
    h : float
        Heat transfer coefficient [W/m^2/K]

    Notes
    -----
    The liquid-only heat transfer coefficient will be calculated with the
    `laminar_entry_Seider_Tate` correlation, should it not be provided as an
    input. Many of the arguments to this function are optional and are only
    used if `hl` is not provided.

    `hl` should be calculated with a velocity equal to that determined with
    a combined volumetric flow of both the liquid and the gas. All other
    parameters used in calculating the heat transfer coefficient are those
    of the liquid. If the viscosity at the wall temperature is not given, the
    liquid viscosity correction in `laminar_entry_Seider_Tate` is not applied.

    Examples
    --------
    >>> Knott(m=1, x=.9, D=.3, rhol=1000, rhog=2.5, Cpl=2300, kl=.6, mu_b=1E-3,
    ... mu_w=1.2E-3, L=4)
    4225.536758045839

    References
    ----------
    .. [1] Knott, R. F., R. N. Anderson, Andreas. Acrivos, and E. E. Petersen.
       "An Experimental Study of Heat Transfer to Nitrogen-Oil Mixtures."
       Industrial & Engineering Chemistry 51, no. 11 (November 1, 1959):
       1369-72. doi:10.1021/ie50599a032.
    .. [2] Dongwoo Kim, Venkata K. Ryali, Afshin J. Ghajar, Ronald L.
       Dougherty. "Comparison of 20 Two-Phase Heat Transfer Correlations with
       Seven Sets of Experimental Data, Including Flow Pattern and Tube
       Inclination Effects." Heat Transfer Engineering 20, no. 1 (February 1,
       1999): 15-40. doi:10.1080/014576399271691.
    '''
    Vgs = m*x/(rhog*pi/4*D**2)
    Vls = m*(1-x)/(rhol*pi/4*D**2)
    if not hl:
        V = Vgs + Vls # Net velocity
        Re = Reynolds(V=V, D=D, rho=rhol, mu=mu_b)
        Pr = Prandtl(Cp=Cpl, k=kl, mu=mu_b)
        Nul = laminar_entry_Seider_Tate(Re=Re, Pr=Pr, L=L, Di=D, mu=mu_b, mu_w=mu_w)
        hl = Nul*kl/D
    return hl*(1 + Vgs/Vls)**(1/3.)
Пример #12
0
def Shah(m, x, D, rhol, mul, kl, Cpl, P, Pc):
    r'''Calculates heat transfer coefficient for condensation
    of a fluid inside a tube, as presented in [1]_ and again by the same
    author in [2]_; also given in [3]_. Requires no properties of the gas.
    Uses the Dittus-Boelter correlation for single phase heat transfer
    coefficient, with a Reynolds number assuming all the flow is liquid.

    .. math::
        h_{TP} = h_L\left[(1-x)^{0.8} +\frac{3.8x^{0.76}(1-x)^{0.04}}
        {P_r^{0.38}}\right]

    Parameters
    ----------
    m : float
        Mass flow rate [kg/s]
    x : float
        Quality at the specific interval [-]
    D : float
        Diameter of the channel [m]
    rhol : float
        Density of the liquid [kg/m^3]
    mul : float
        Viscosity of liquid [Pa*s]
    kl : float
        Thermal conductivity of liquid [W/m/K]
    Cpl : float
        Constant-pressure heat capacity of liquid [J/kg/K]
    P : float
        Pressure of the fluid, [Pa]
    Pc : float
        Critical pressure of the fluid, [Pa]

    Returns
    -------
    h : float
        Heat transfer coefficient [W/m^2/K]

    Notes
    -----
    [1]_ is well written an unambiguous as to how to apply this equation.

    Examples
    --------
    >>> Shah(m=1, x=0.4, D=.3, rhol=800, mul=1E-5, kl=0.6, Cpl=2300, P=1E6, Pc=2E7)
    2561.2593415479214

    References
    ----------
    .. [1] Shah, M. M. "A General Correlation for Heat Transfer during Film
       Condensation inside Pipes." International Journal of Heat and Mass
       Transfer 22, no. 4 (April 1, 1979): 547-56.
       doi:10.1016/0017-9310(79)90058-9.
    .. [2] Shah, M. M., Heat Transfer During Film Condensation in Tubes and
       Annuli: A Review of the Literature, ASHRAE Transactions, vol. 87, no.
       3, pp. 1086-1100, 1981.
    .. [3] Kakaç, Sadik, ed. Boilers, Evaporators, and Condensers. 1st.
       Wiley-Interscience, 1991.
    '''
    VL = m / (rhol * pi / 4 * D**2)
    ReL = Reynolds(V=VL, D=D, rho=rhol, mu=mul)
    Prl = Prandtl(Cp=Cpl, k=kl, mu=mul)
    hL = turbulent_Dittus_Boelter(ReL, Prl) * kl / D
    Pr = P / Pc
    return hL * ((1 - x)**0.8 + 3.8 * x**0.76 * (1 - x)**0.04 / Pr**0.38)
Пример #13
0
def Cavallini_Smith_Zecchin(m, x, D, rhol, rhog, mul, mug, kl, Cpl):
    r'''Calculates heat transfer coefficient for condensation
    of a fluid inside a tube, as presented in
    [1]_, also given in [2]_ and [3]_.

    .. math::
        Nu = \frac{hD_i}{k_l} = 0.05 Re_e^{0.8} Pr_l^{0.33}

    .. math::
        Re_{eq} = Re_g(\mu_g/\mu_l)(\rho_l/\rho_g)^{0.5} + Re_l

    .. math::
        v_{gs} = \frac{mx}{\rho_g \frac{\pi}{4}D^2}

    .. math::
        v_{ls} = \frac{m(1-x)}{\rho_l \frac{\pi}{4}D^2}

    Parameters
    ----------
    m : float
        Mass flow rate [kg/s]
    x : float
        Quality at the specific interval [-]
    D : float
        Diameter of the channel [m]
    rhol : float
        Density of the liquid [kg/m^3]
    rhog : float
        Density of the gas [kg/m^3]
    mul : float
        Viscosity of liquid [Pa*s]
    mug : float
        Viscosity of gas [Pa*s]
    kl : float
        Thermal conductivity of liquid [W/m/K]
    Cpl : float
        Constant-pressure heat capacity of liquid [J/kg/K]

    Returns
    -------
    h : float
        Heat transfer coefficient [W/m^2/K]

    Notes
    -----

    Examples
    --------
    >>> Cavallini_Smith_Zecchin(m=1, x=0.4, D=.3, rhol=800, rhog=2.5, mul=1E-5, mug=1E-3, kl=0.6, Cpl=2300)
    5578.218369177804

    References
    ----------
    .. [1] A. Cavallini, J. R. Smith and R. Zecchin, A dimensionless correlation
       for heat transfer in forced convection condensation, 6th International
       Heat Transfer Conference., Tokyo, Japan (1974) 309-313.
    .. [2] Kakaç, Sadik, ed. Boilers, Evaporators, and Condensers. 1st.
       Wiley-Interscience, 1991.
    .. [3] Balcılar, Muhammet, Ahmet Selim Dalkılıç, Berna Bolat, and Somchai
       Wongwises. "Investigation of Empirical Correlations on the Determination
       of Condensation Heat Transfer Characteristics during Downward Annular
       Flow of R134a inside a Vertical Smooth Tube Using Artificial
       Intelligence Algorithms." Journal of Mechanical Science and Technology
       25, no. 10 (October 12, 2011): 2683-2701. doi:10.1007/s12206-011-0618-2.
    '''
    Prl = Prandtl(Cp=Cpl, mu=mul, k=kl)
    Vl = m * (1 - x) / (rhol * pi / 4 * D**2)
    Vg = m * x / (rhog * pi / 4 * D**2)
    Rel = Reynolds(V=Vl, D=D, rho=rhol, mu=mul)
    Reg = Reynolds(V=Vg, D=D, rho=rhog, mu=mug)
    '''The following was coded, and may be used instead of the above lines,
    to check that the definitions of parameters here provide the same results
    as those defined in [1]_.
    G = m/(pi/4*D**2)
    Re = G*D/mul
    Rel = Re*(1-x)
    Reg = Re*x/(mug/mul)'''
    Reeq = Reg * (mug / mul) * (rhol / rhog)**0.5 + Rel
    Nul = 0.05 * Reeq**0.8 * Prl**0.33
    return Nul * kl / D  # confirmed to be with respect to the liquid