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.)
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.)
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)
def accept(self): Dp=Ltot=nc=0 elements=list() Q=float(self.form.editFlow.text())/3600 if not self.isLiquid: Q=Q/self.Rho if self.form.comboWhat.currentText()=='<on selection>': elements = FreeCADGui.Selection.getSelection() else: o=FreeCAD.ActiveDocument.getObjectsByLabel(self.form.comboWhat.currentText())[0] if hasattr(o,'PType') and o.PType=='PypeBranch': elements=[FreeCAD.ActiveDocument.getObject(name) for name in o.Tubes+o.Curves] elif hasattr(o,'PType') and o.PType=='PypeLine': group=FreeCAD.ActiveDocument.getObjectsByLabel(o.Label+'_pieces')[0] elements=group.OutList self.form.editResults.clear() for o in elements: loss=0 if hasattr(o,'PType') and o.PType in ['Pipe','Elbow']: ID=float(o.ID)/1000 e=float(self.form.editRough.text())*1e-6/ID v=Q/((ID)**2*pi/4) Re=Reynolds(V=v,D=ID,rho=self.Rho, mu=self.Mu) f=friction.friction_factor(Re, eD=e) # Darcy, =4xFanning if o.PType=='Pipe': L=float(o.Height)/1000 Ltot+=L K=K_from_f(fd=f, L=L, D=ID) loss=dP_from_K(K,rho=self.Rho,V=v) FreeCAD.Console.PrintMessage('%s: %s\nID=%.2f\nV=%.2f\ne=%f\nf=%f\nK=%f\nL=%.3f\nDp = %.5f bar\n***'%(o.PType,o.Label,ID,v,e,f,K,L,loss/1e5)) self.form.editResults.append('%s\t%.1f mm\t%.1f m/s\t%.5f bar'%(o.Label,ID*1000,v,loss/1e5)) elif o.PType=='Elbow': ang=float(o.BendAngle) R=float(o.BendRadius)/1000 nc+=1 K=fittings.bend_rounded(ID,ang,f,R) loss=dP_from_K(K,rho=self.Rho,V=v) FreeCAD.Console.PrintMessage('%s: %s\nID=%.2f\nV=%.2f\ne=%f\nf=%f\nK=%f\nang=%.3f\nR=%f\nDp = %.5f bar\n***'%(o.PType,o.Label,ID,v,e,f,K,ang,R,loss/1e5)) self.form.editResults.append('%s\t%.1f mm\t%.1f m/s\t%.5f bar'%(o.Label,ID*1000,v,loss/1e5)) elif o.PType=='Reduct': pass elif hasattr(o,'Kv') and o.Kv>0: if self.isLiquid: loss=(Q*3600/o.Kv)**2*100000 elif self.form.comboFluid.currentText()=='water' and not self.isLiquid: pass # TODO formula for steam else: pass # TODO formula for gases if hasattr(o,'ID'): ID = float(o.ID)/1000 v=Q/(ID**2*pi/4) else: v = 0 ID = 0 self.form.editResults.append('%s\t%.1f mm\t%.1f m/s\t%.5f bar'%(o.Label,ID*1000,v,loss/1e5)) Dp+=loss if Dp>200: result=' = %.3f bar'%(Dp/100000) else: result=' = %.2e bar'%(Dp/100000) self.form.labResult.setText(result) self.form.labLength.setText('Total length = %.3f m' %Ltot) self.form.labCurves.setText('Nr. of curves = %i' %nc)
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)
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} Re_{eq} = Re_g(\mu_g/\mu_l)(\rho_l/\rho_g)^{0.5} + Re_l v_{gs} = \frac{mx}{\rho_g \frac{\pi}{4}D^2} 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