def Prandtl_von_Karman_Nikuradse_numeric(Re): rat = 2.51 / Re def to_solve(f): # Good to 1E75, down to 1E-17 v = f**-0.5 return v + 2.0 * log10(rat * v) return secant(to_solve, 0.000001)
def test_valve_coefficients(): Cv = Kv_to_Cv(2) assert_close(Cv, 2.3121984567073133) Kv = Cv_to_Kv(2.312) assert_close(Kv, 1.9998283393826013) K = Kv_to_K(2.312, .015) assert_close(K, 15.15337460039990) Kv = K_to_Kv(15.15337460039990, .015) assert_close(Kv, 2.312) # Two way conversions K = Cv_to_K(2.712, .015) assert_close(K, 14.719595348352552) assert_close(K, Kv_to_K(Cv_to_Kv(2.712), 0.015)) Cv = K_to_Cv(14.719595348352552, .015) assert_close(Cv, 2.712) assert_close(Cv, Kv_to_Cv(K_to_Kv(14.719595348352552, 0.015))) # Code to generate the Kv Cv conversion factor # Round 1 trip; randomly assume Kv = 12, rho = 900; they can be anything # an tit still works dP = 1E5 rho = 900. Kv = 12. Q = Kv/3600. D = .01 V = Q/(pi/4*D**2) K = dP/(.5*rho*V*V) good_K = K def to_solve(x): from fluids.constants import gallon, minute, hour, psi conversion = gallon/minute*hour # from gpm to m^3/hr dP = 1*psi Cv = Kv*x*conversion Q = Cv/3600 D = .01 V = Q/(pi/4*D**2) K = dP/(.5*rho*V*V) return K - good_K ans = secant(to_solve, 1.2) assert_close(ans, 1.1560992283536566)
def test_cone_meter_expansibility_Stewart_full(): err = lambda Dc, beta : diameter_ratio_cone_meter(D=1, Dc=Dc) - beta solve_Dc = lambda beta : float(secant(err, .7, args=(beta,))) # Accidentally missed the beta ratio 0.75, oops vals = [[1.0000, 0.9887, 0.9774, 0.9661, 0.9548, 0.9435, 0.9153, 0.8871, 0.8588], [1.0000, 0.9885, 0.9769, 0.9654, 0.9538, 0.9423, 0.9134, 0.8846, 0.8557], [1.0000, 0.9881, 0.9762, 0.9644, 0.9525, 0.9406, 0.9109, 0.8812, 0.8515], [1.0000, 0.9877, 0.9754, 0.9630, 0.9507, 0.9384, 0.9076, 0.8768, 0.8460], [1.0000, 0.9871, 0.9742, 0.9613, 0.9485, 0.9356, 0.9033, 0.8711, 0.8389], [1.0000, 0.9864, 0.9728, 0.9592, 0.9456, 0.9320, 0.8980, 0.8640, 0.8300]] pressure_ratios = [1, 0.98, 0.96, 0.94, 0.92, 0.9, 0.85, 0.8, 0.75] betas = [.45, .5, .55, .6, .65, .7, .75] k = 1.2 for i, beta in enumerate(betas[:-1]): Dc = solve_Dc(beta) for j, pr in enumerate(pressure_ratios): eps = cone_meter_expansibility_Stewart(D=1, Dc=Dc, P1=1E5, P2=pr*1E5, k=1.2) eps = round(eps, 4) assert eps == vals[i][j]
def Thome(m, x, D, rhol, rhog, mul, mug, kl, kg, Cpl, Cpg, Hvap, sigma, Psat, Pc, q=None, Te=None): r'''Calculates heat transfer coefficient for film boiling of saturated fluid in any orientation of flow. Correlation is as developed in [1]_ and [2]_, and also reviewed [3]_. This is a complicated model, but expected to have more accuracy as a result. Either the heat flux or excess temperature is required for the calculation of heat transfer coefficient. The solution for a specified excess temperature is solved numerically, making it slow. .. math:: h(z) = \frac{t_l}{\tau} h_l(z) +\frac{t_{film}}{\tau} h_{film}(z) + \frac{t_{dry}}{\tau} h_{g}(z) .. math:: h_{l/g}(z) = (Nu_{lam}^4 + Nu_{trans}^4)^{1/4} k/D .. math:: Nu_{laminar} = 0.91 {Pr}^{1/3} \sqrt{ReD/L(z)} .. math:: Nu_{trans} = \frac{ (f/8) (Re-1000)Pr}{1+12.7 (f/8)^{1/2} (Pr^{2/3}-1)} \left[ 1 + \left( \frac{D}{L(z)}\right)^{2/3}\right] .. math:: f = (1.82 \log_{10} Re - 1.64 )^{-2} .. math:: L_l = \frac{\tau G_{tp}}{\rho_l}(1-x) .. math:: L_{dry} = v_p t_{dry} .. math:: t_l = \frac{\tau}{1 + \frac{\rho_l}{\rho_g}\frac{x}{1-x}} .. math:: t_v = \frac{\tau}{1 + \frac{\rho_g}{\rho_l}\frac{1-x}{x}} .. math:: \tau = \frac{1}{f_{opt}} .. math:: f_{opt} = \left(\frac{q}{q_{ref}}\right)^{n_f} .. math:: q_{ref} = 3328\left(\frac{P_{sat}}{P_c}\right)^{-0.5} .. math:: t_{dry,film} = \frac{\rho_l \Delta H_{vap}}{q}[\delta_0(z) - \delta_{min}] .. math:: \frac{\delta_0}{D} = C_{\delta 0}\left(3\sqrt{\frac{\nu_l}{v_p D}} \right)^{0.84}\left[(0.07Bo^{0.41})^{-8} + 0.1^{-8}\right]^{-1/8} .. math:: Bo = \frac{\rho_l D}{\sigma} v_p^2 .. math:: v_p = G_{tp} \left[\frac{x}{\rho_g} + \frac{1-x}{\rho_l}\right] .. math:: h_{film}(z) = \frac{2 k_l}{\delta_0(z) + \delta_{min}(z)} .. math:: \delta_{min} = 0.3\cdot 10^{-6} \text{m} .. math:: C_{\delta,0} = 0.29 .. math:: n_f = 1.74 if t dry film > tv: .. math:: \delta_{end}(x) = \delta(z, t_v) .. math:: t_{film} = t_v .. math:: t_{dry} = 0 Otherwise: .. math:: \delta_{end}(z) = \delta_{min} .. math:: t_{film} = t_{dry,film} .. math:: t_{dry} = t_v - t_{film} 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] mul : float Viscosity of liquid [Pa*s] mug : float Viscosity of gas [Pa*s] kl : float Thermal conductivity of liquid [W/m/K] kg : float Thermal conductivity of gas [W/m/K] Cpl : float Heat capacity of liquid [J/kg/K] Cpg : float Heat capacity of gas [J/kg/K] Hvap : float Heat of vaporization of liquid [J/kg] sigma : float Surface tension of liquid [N/m] Psat : float Vapor pressure of fluid, [Pa] Pc : float Critical pressure of fluid, [Pa] q : float, optional Heat flux to wall [W/m^2] Te : float, optional Excess temperature of wall, [K] Returns ------- h : float Heat transfer coefficient [W/m^2/K] Notes ----- [1]_ and [2]_ have been reviewed, and are accurately reproduced in [3]_. [1]_ used data from 7 studies, covering 7 fluids and Dh from 0.7-3.1 mm, heat flux from 0.5-17.8 W/cm^2, x from 0.01-0.99, and G from 50-564 kg/m^2/s. Liquid and/or gas slugs are both considered, and are hydrodynamically developing. `Ll` is the calculated length of liquid slugs, and `L_dry` is the same for vapor slugs. Because of the complexity of the model and that there is some logic in this function, `Te` as an input may lead to a different solution that the calculated `q` will in return. Examples -------- >>> Thome(m=1, x=0.4, D=0.3, rhol=567., rhog=18.09, kl=0.086, kg=0.2, ... mul=156E-6, mug=1E-5, Cpl=2300, Cpg=1400, sigma=0.02, Hvap=9E5, ... Psat=1E5, Pc=22E6, q=1E5) 1633.008836502032 References ---------- .. [1] Thome, J. R., V. Dupont, and A. M. Jacobi. "Heat Transfer Model for Evaporation in Microchannels. Part I: Presentation of the Model." International Journal of Heat and Mass Transfer 47, no. 14-16 (July 2004): 3375-85. doi:10.1016/j.ijheatmasstransfer.2004.01.006. .. [2] Dupont, V., J. R. Thome, and A. M. Jacobi. "Heat Transfer Model for Evaporation in Microchannels. Part II: Comparison with the Database." International Journal of Heat and Mass Transfer 47, no. 14-16 (July 2004): 3387-3401. doi:10.1016/j.ijheatmasstransfer.2004.01.007. .. [3] Bertsch, Stefan S., Eckhard A. Groll, and Suresh V. Garimella. "Review and Comparative Analysis of Studies on Saturated Flow Boiling in Small Channels." Nanoscale and Microscale Thermophysical Engineering 12, no. 3 (September 4, 2008): 187-227. doi:10.1080/15567260802317357. ''' if q is None and Te is not None: q = secant(to_solve_q_Thome, 1E4, args=( m, x, D, rhol, rhog, kl, kg, mul, mug, Cpl, Cpg, sigma, Hvap, Psat, Pc, Te)) return Thome(m=m, x=x, D=D, rhol=rhol, rhog=rhog, kl=kl, kg=kg, mul=mul, mug=mug, Cpl=Cpl, Cpg=Cpg, sigma=sigma, Hvap=Hvap, Psat=Psat, Pc=Pc, q=q) elif q is None and Te is None: raise ValueError('Either q or Te is needed for this correlation') C_delta0 = 0.3E-6 G = m/(pi/4*D**2) Rel = G*D*(1-x)/mul Reg = G*D*x/mug qref = 3328*(Psat/Pc)**-0.5 if q is None: q = 1e4 # Make numba happy, their bug, never gets ran fopt = (q/qref)**1.74 tau = 1./fopt vp = G*(x/rhog + (1-x)/rhol) Bo = rhol*D/sigma*vp**2 # Not standard definition nul = mul/rhol delta0 = D*0.29*(3*(nul/vp/D)**0.5)**0.84*((0.07*Bo**0.41)**-8 + 0.1**-8)**(-1/8.) tl = tau/(1 + rhol/rhog*(x/(1.-x))) tv = tau/(1 ++ rhog/rhol*((1.-x)/x)) t_dry_film = rhol*Hvap/q*(delta0 - C_delta0) if t_dry_film > tv: t_film = tv delta_end = delta0 - q/rhol/Hvap*tv # what could time possibly be? t_dry = 0 else: t_film = t_dry_film delta_end = C_delta0 t_dry = tv-t_film Ll = tau*G/rhol*(1-x) Ldry = t_dry*vp Prg = Prandtl(Cp=Cpg, k=kg, mu=mug) Prl = Prandtl(Cp=Cpl, k=kl, mu=mul) fg = (1.82*log10(Reg) - 1.64)**-2 fl = (1.82*log10(Rel) - 1.64)**-2 Nu_lam_Zl = 2*0.455*(Prl)**(1/3.)*(D*Rel/Ll)**0.5 Nu_trans_Zl = turbulent_Gnielinski(Re=Rel, Pr=Prl, fd=fl)*(1 + (D/Ll)**(2/3.)) if Ldry == 0: Nu_lam_Zg, Nu_trans_Zg = 0, 0 else: Nu_lam_Zg = 2*0.455*(Prg)**(1/3.)*(D*Reg/Ldry)**0.5 Nu_trans_Zg = turbulent_Gnielinski(Re=Reg, Pr=Prg, fd=fg)*(1 + (D/Ldry)**(2/3.)) h_Zg = kg/D*(Nu_lam_Zg**4 + Nu_trans_Zg**4)**0.25 h_Zl = kl/D*(Nu_lam_Zl**4 + Nu_trans_Zl**4)**0.25 h_film = 2*kl/(delta0 + C_delta0) return tl/tau*h_Zl + t_film/tau*h_film + t_dry/tau*h_Zg
def Nu_Nusselt_Rayleigh_Holling_Herwig(Pr, Gr, buoyancy=True): r'''Calculates the Nusselt number for natural convection between two theoretical flat plates. The height between the plates is infinite, and one of the other dimensions of the plates is much larger than the other. This model is a non-linear equation which is solved numerically. The model can calculate `Nu` for `Ra` ranges between 350 and larger numbers; [1]_ recommends :math:`10^{5} < Ra < 10^{15}`. .. math:: \text{Nu} = \frac{{Ra}^{1/3}}{[0.05\log(\frac{0.078}{16}{Ra}^{1.323}) + 2D]^{4/3}} .. math:: D = -\frac{14.94}{{Ra}^{0.25}} + 3.43 Parameters ---------- Pr : float Prandtl number with respect to fluid properties [-] Gr : float Grashof number with respect to fluid properties and plate - plate temperature difference [-] buoyancy : bool, optional Whether or not the plate's free convection is buoyancy assisted (hot plate) or not, [-] Returns ------- Nu : float Nusselt number with respect to height between the two plates, [-] Notes ----- A range of calculated values are provided in [1]_; they all match the results of this function. This model is recommended in [2]_. For :math:`Ra < 1708`, `Nu` = 1; for cases not assited by `buoyancy`, `Nu` is also 1. Examples -------- >>> Nu_Nusselt_Rayleigh_Holling_Herwig(5.54, 3.21e8, buoyancy=True) 77.54656801896913 References ---------- .. [1] Hölling, M., and H. Herwig. "Asymptotic Analysis of Heat Transfer in Turbulent Rayleigh–Bénard Convection." International Journal of Heat and Mass Transfer 49, no. 5 (March 1, 2006): 1129-36. https://doi.org/10.1016/j.ijheatmasstransfer.2005.09.002. .. [2] Gesellschaft, V. D. I., ed. VDI Heat Atlas. 2nd ed. 2010 edition. Berlin ; New York: Springer, 2010. ''' if not buoyancy: return 1.0 Rac = 1708 # Constant Ra = Gr * Pr if Ra < Rac: return 1.0 Ra_third = Ra**(1.0 / 3.0) D2 = 2.0 * (-14.94 * Ra**-0.25 + 3.43) def to_solve(Nu): err = Ra_third * (0.1 / 2.0 * log(1.0 / 16.0 * Ra * Nu) + D2)**(-4.0 / 3.0) - Nu return err Nu_guess = Ra_third * (0.1 / 2.0 * log(.078 / 16.0 * Ra**1.323) + D2)**(-4.0 / 3.0) return secant(to_solve, Nu_guess)
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))
def Stichlmair_wet(Vg, Vl, rhog, rhol, mug, voidage, specific_area, C1, C2, C3, H=1): r'''Calculates dry pressure drop across a packed column, using the Stichlmair [1]_ correlation. Uses three regressed constants for each type of packing, and voidage and specific area. This model is for irrigated columns only. Pressure drop is given by: .. math:: \frac{\Delta P_{irr}}{H} = \frac{\Delta P_{dry}}{H}\left(\frac {1-\epsilon + h_T}{1-\epsilon}\right)^{(2+c)/3} \left(\frac{\epsilon}{\epsilon-h_T}\right)^{4.65} .. math:: h_T = h_0\left[1 + 20\left(\frac{\Delta Pirr}{H\rho_L g}\right)^2\right] .. math:: Fr_L = \frac{V_L^2 a}{g \epsilon^{4.65}} .. math:: h_0 = 0.555 Fr_L^{1/3} .. math:: c = \frac{-C_1/Re_g - C_2/(2Re_g^{0.5})}{f_0} .. math:: \Delta P_{dry} = \frac{3}{4} f_0 \frac{1-\epsilon}{\epsilon^{4.65}} \rho_G \frac{H}{d_p}V_g^2 .. math:: f_0 = \frac{C_1}{Re_g} + \frac{C_2}{Re_g^{0.5}} + C_3 .. math:: d_p = \frac{6(1-\epsilon)}{a} Parameters ---------- Vg : float Superficial velocity of gas, Q/A [m/s] Vl : float Superficial velocity of liquid, Q/A [m/s] rhog : float Density of gas [kg/m^3] rhol : float Density of liquid [kg/m^3] mug : float Viscosity of gas [Pa*s] voidage : float Voidage of bed of packing material [] specific_area : float Specific area of the packing material [m^2/m^3] C1 : float Packing-specific constant [] C2 : float Packing-specific constant [] C3 : float Packing-specific constant [] H : float, optional Height of packing [m] Returns ------- dP : float Pressure drop across irrigated packing [Pa] Notes ----- This model is used by most process simulation tools. If H is not provided, it defaults to 1. If Z is not provided, it defaults to 1. A numerical solver is used and needed by this model. Its initial guess is the dry pressure drop. Convergence problems may occur. The model as described in [1]_ appears to have a typo, and could not match the example. As described in [2]_, however, the model works. Examples -------- Example is from [1]_, matches. >>> Stichlmair_wet(Vg=0.4, Vl = 5E-3, rhog=5., rhol=1200., mug=5E-5, ... voidage=0.68, specific_area=260., C1=32., C2=7., C3=1.) 539.8768237253518 References ---------- .. [1] Stichlmair, J., J. L. Bravo, and J. R. Fair. "General Model for Prediction of Pressure Drop and Capacity of Countercurrent Gas/liquid Packed Columns." Gas Separation & Purification 3, no. 1 (March 1989): 19-28. doi:10.1016/0950-4214(89)80016-7. .. [2] Piche, Simon R., Faical Larachi, and Bernard P. A. Grandjean. "Improving the Prediction of Irrigated Pressure Drop in Packed Absorption Towers." The Canadian Journal of Chemical Engineering 79, no. 4 (August 1, 2001): 584-94. doi:10.1002/cjce.5450790417. ''' dp = 6.0 * (1.0 - voidage) / specific_area Re = Vg * rhog * dp / mug f0 = C1 / Re + C2 / Re**0.5 + C3 dP_dry = 3 / 4. * f0 * (1 - voidage) / voidage**4.65 * rhog * H / dp * Vg * Vg c = (-C1 / Re - C2 / (2 * Re**0.5)) / f0 Frl = Vl**2 * specific_area / (g * voidage**4.65) h0 = 0.555 * Frl**(1 / 3.) c1 = 1.0 / (H * rhol * g) c1 *= c1 return secant(_Stichlmair_wet_err, dP_dry, args=(h0, c1, dP_dry, H, voidage, c))
def liquid_jet_pump_ancillary(rhop, rhos, Kp, Ks, d_nozzle=None, d_mixing=None, Qp=None, Qs=None, P1=None, P2=None): r'''Calculates the remaining variable in a liquid jet pump when solving for one if the inlet variables only and the rest of them are known. The equation comes from conservation of energy and momentum in the mixing chamber. The variable to be solved for must be one of `d_nozzle`, `d_mixing`, `Qp`, `Qs`, `P1`, or `P2`. .. math:: P_1 - P_2 = \frac{1}{2}\rho_pV_n^2(1+K_p) - \frac{1}{2}\rho_s V_3^2(1+K_s) Rearrange to express V3 in terms of Vn, and using the density ratio `C`, the expression becomes: .. math:: P_1 - P_2 = \frac{1}{2}\rho_p V_n^2\left[(1+K_p) - C(1+K_s) \left(\frac{MR}{1-R}\right)^2\right] Using the primary nozzle area and flow rate: .. math:: P_1 - P_2 = \frac{1}{2}\rho_p \left(\frac{Q_p}{A_n}\right)^2 \left[(1+K_p) - C(1+K_s) \left(\frac{MR}{1-R}\right)^2\right] For `P`, `P2`, `Qs`, and `Qp`, the equation can be rearranged explicitly for them. For `d_mixing` and `d_nozzle`, a bounded solver is used searching between 1E-9 m and 20 times the other diameter which was specified. Parameters ---------- rhop : float The density of the primary (motive) fluid, [kg/m^3] rhos : float The density of the secondary fluid (drawn from the vacuum chamber), [kg/m^3] Kp : float The primary nozzle loss coefficient, [-] Ks : float The secondary inlet loss coefficient, [-] d_nozzle : float, optional The inside diameter of the primary fluid's nozle, [m] d_mixing : float, optional The diameter of the mixing chamber, [m] Qp : float, optional The volumetric flow rate of the primary fluid, [m^3/s] Qs : float, optional The volumetric flow rate of the secondary fluid, [m^3/s] P1 : float, optional The pressure of the primary fluid entering its nozzle, [Pa] P2 : float, optional The pressure of the secondary fluid at the entry of the ejector, [Pa] Returns ------- solution : float The parameter not specified (one of `d_nozzle`, `d_mixing`, `Qp`, `Qs`, `P1`, or `P2`), (units of `m`, `m`, `m^3/s`, `m^3/s`, `Pa`, or `Pa` respectively) Notes ----- The following SymPy code was used to obtain the analytical formulas ( they are not shown here due to their length): >>> from sympy import * >>> A_nozzle, A_mixing, Qs, Qp, P1, P2, rhos, rhop, Ks, Kp = symbols('A_nozzle, A_mixing, Qs, Qp, P1, P2, rhos, rhop, Ks, Kp') >>> R = A_nozzle/A_mixing >>> M = Qs/Qp >>> C = rhos/rhop >>> rhs = rhop/2*(Qp/A_nozzle)**2*((1+Kp) - C*(1 + Ks)*((M*R)/(1-R))**2 ) >>> new = Eq(P1 - P2, rhs) >>> #solve(new, Qp) >>> #solve(new, Qs) >>> #solve(new, P1) >>> #solve(new, P2) Examples -------- Calculating primary fluid nozzle inlet pressure P1: >>> liquid_jet_pump_ancillary(rhop=998., rhos=1098., Ks=0.11, Kp=.04, ... P2=133600, Qp=0.01, Qs=0.01, d_mixing=0.045, d_nozzle=0.02238) 426434.60314398084 References ---------- .. [1] Ejectors and Jet Pumps. Design and Performance for Incompressible Liquid Flow. 85032. ESDU International PLC, 1985. ''' unknowns = sum(i is None for i in (d_nozzle, d_mixing, Qs, Qp, P1, P2)) if unknowns > 1: raise ValueError('Too many unknowns') elif unknowns < 1: raise ValueError('Overspecified') C = rhos / rhop if Qp is not None and Qs is not None: M = Qs / Qp if d_nozzle is not None: A_nozzle = pi / 4 * d_nozzle * d_nozzle if d_mixing is not None: A_mixing = pi / 4 * d_mixing * d_mixing R = A_nozzle / A_mixing if P1 is None: return rhop / 2 * (Qp / A_nozzle)**2 * ((1 + Kp) - C * (1 + Ks) * ((M * R) / (1 - R))**2) + P2 elif P2 is None: return -rhop / 2 * (Qp / A_nozzle)**2 * ((1 + Kp) - C * (1 + Ks) * ((M * R) / (1 - R))**2) + P1 elif Qs is None: try: return sqrt( (-2 * A_nozzle**2 * P1 + 2 * A_nozzle**2 * P2 + Kp * Qp**2 * rhop + Qp**2 * rhop) / (C * rhop * (Ks + 1))) * (A_mixing - A_nozzle) / A_nozzle except ValueError: return -1j elif Qp is None: return A_nozzle * sqrt( (2 * A_mixing**2 * P1 - 2 * A_mixing**2 * P2 - 4 * A_mixing * A_nozzle * P1 + 4 * A_mixing * A_nozzle * P2 + 2 * A_nozzle**2 * P1 - 2 * A_nozzle**2 * P2 + C * Ks * Qs**2 * rhop + C * Qs**2 * rhop) / (rhop * (Kp + 1))) / (A_mixing - A_nozzle) elif d_nozzle is None: def err(d_nozzle): return P1 - liquid_jet_pump_ancillary(rhop=rhop, rhos=rhos, Kp=Kp, Ks=Ks, d_nozzle=d_nozzle, d_mixing=d_mixing, Qp=Qp, Qs=Qs, P1=None, P2=P2) return brenth(err, 1E-9, d_mixing * 20) elif d_mixing is None: def err(d_mixing): return P1 - liquid_jet_pump_ancillary(rhop=rhop, rhos=rhos, Kp=Kp, Ks=Ks, d_nozzle=d_nozzle, d_mixing=d_mixing, Qp=Qp, Qs=Qs, P1=None, P2=P2) try: return brenth(err, 1E-9, d_nozzle * 20) except: return secant(err, d_nozzle * 2)
def solve_PT_HSGUA_NP_guess_bisect(self, zs, fixed_val, spec_val, fixed_var='P', spec='H', iter_var='T'): phases = self.phases constants = self.constants correlations = self.correlations min_bound, max_bound = self.bounds_PT_HSGUA() init_methods = [SHAW_ELEMENTAL, IDEAL_WILSON] guess = None for method in init_methods: try: guess, VF, xs, ys = TPV_solve_HSGUA_guesses_VL( zs, method, constants, correlations, fixed_val, spec_val, iter_var=iter_var, fixed_var=fixed_var, spec=spec, maxiter=50, xtol=1E-5, ytol=None, bounded=False, min_bound=min_bound, max_bound=max_bound, user_guess=None, last_conv=None, T_ref=298.15, P_ref=101325.0) break except NotImplementedError: continue except Exception as e: #print(e) pass if guess is None: if iter_var == 'T': guess = 298.15 elif iter_var == 'P': guess = 101325.0 elif iter_var == 'V': guess = 0.024465403697038125 sln = [] global iterations iterations = 0 kwargs = {fixed_var: fixed_val, 'zs': zs} def to_solve(iter_val): global iterations iterations += 1 kwargs[iter_var] = iter_val res = self.flash(**kwargs) err = getattr(res, spec)() - spec_val sln[:] = (res, iter_val) return err ytol = abs(spec_val) * self.TPV_HSGUA_BISECT_YTOL sln_val = secant(to_solve, guess, xtol=self.TPV_HSGUA_BISECT_XTOL, ytol=ytol, require_xtol=self.TPV_HSGUA_BISECT_YTOL_ONLY, require_eval=True, bisection=True, low=min_bound, high=max_bound) return sln[0], {'iterations': iterations, 'err': sln[1]}
def flash_ideal(zs, funcs, Tcs=None, T=None, P=None, VF=None): r'''PVT flash model using ideal, composition-independent equation. Solves the various cases of composition-independent models. Capable of solving with two of `T`, `P`, and `VF` for the other one; that results in three solve modes, but for `VF=1` and `VF=0`, there are additional solvers; for a total of seven solvers implemented. The function takes a list of callables that take `T` in Kelvin as an argument, and return vapor pressure. The callables can include the effect of non-ideal pure component fugacity coefficients. For the (`T`, `P`) and (`P`, `VF`) cases, the Poynting correction factor can be easily included as well but not the (`T`, `VF`) case as the callable only takes `T` as an argument. Normally the Poynting correction factor is used with activity coefficient models with composition dependence. Both `flash_wilson` and `flash_Tb_Tc_Pc` are specialized cases of this function and have the same functionality but with the model built right in. Even when using more complicated models, this is useful for obtaining initial This model uses `flash_inner_loop` to solve the Rachford-Rice problem. Parameters ---------- zs : list[float] Mole fractions of the phase being flashed, [-] funcs : list[Callable] Functions to calculate ideal or real vapor pressures, take temperature in Kelvin and return pressure in Pa, [-] Tcs : list[float], optional Critical temperatures of all species; uses as upper bounds and only for the case that `T` is not specified; if they are needed and not given, it is assumed a method `solve_prop` exists in each of `funcs` which will accept `P` in Pa and return temperature in `K`, [K] T : float, optional Temperature, [K] P : float, optional Pressure, [Pa] VF : float, optional Molar vapor fraction, [-] Returns ------- T : float Temperature, [K] P : float Pressure, [Pa] VF : float Molar vapor fraction, [-] xs : list[float] Mole fractions of liquid phase, [-] ys : list[float] Mole fractions of vapor phase, [-] Notes ----- For the cases where `VF` is 1 or 0 and T is known, an explicit solution is used. For the same cases where `P` and `VF` are known, there is no explicit solution available. There is an internal `Tmax` parameter, set to 50000 K; which, in the event of convergence of the Secant method, is used as a bounded for a bounded solver. It is used in the PVF solvers. Examples -------- Basic case with four compounds, usingthe Antoine equation as a model and solving for vapor pressure: >>> from chemicals import Antoine, Ambrose_Walton >>> Tcs = [369.83, 425.12, 469.7, 507.6] >>> Antoine_As = [8.92828, 8.93266, 8.97786, 9.00139] >>> Antoine_Bs = [803.997, 935.773, 1064.84, 1170.88] >>> Antoine_Cs = [-26.11, -34.361, -41.136, -48.833] >>> Psat_funcs = [] >>> for i in range(4): ... def Psat_func(T, A=Antoine_As[i], B=Antoine_Bs[i], C=Antoine_Cs[i]): ... return Antoine(T, A, B, C) ... Psat_funcs.append(Psat_func) >>> zs = [.4, .3, .2, .1] >>> T, P, VF, xs, ys = flash_ideal(T=330.55, P=1e6, zs=zs, funcs=Psat_funcs, Tcs=Tcs) >>> round(VF, 10) 1.00817e-05 Similar case, using the Ambrose-Walton corresponding states method to estimate vapor pressures: >>> Tcs = [369.83, 425.12, 469.7, 507.6] >>> Pcs = [4248000.0, 3796000.0, 3370000.0, 3025000.0] >>> omegas = [0.152, 0.193, 0.251, 0.2975] >>> Psat_funcs = [] >>> for i in range(4): ... def Psat_func(T, Tc=Tcs[i], Pc=Pcs[i], omega=omegas[i]): ... return Ambrose_Walton(T, Tc, Pc, omega) ... Psat_funcs.append(Psat_func) >>> _, P, VF, xs, ys = flash_ideal(T=329.151, VF=0, zs=zs, funcs=Psat_funcs, Tcs=Tcs) >>> round(P, 3) 1000013.343 Case with fugacities in the liquid phase, vapor phase, activity coefficients in the liquid phase, and Poynting correction factors. >>> Tcs = [647.14, 514.0] >>> Antoine_As = [10.1156, 10.3368] >>> Antoine_Bs = [1687.54, 1648.22] >>> Antoine_Cs = [-42.98, -42.232] >>> gammas = [1.1, .75] >>> fugacities_gas = [.995, 0.98] >>> fugacities_liq = [.9999, .9998] >>> Poyntings = [1.000001, .999999] >>> zs = [.5, .5] >>> funcs = [] >>> for i in range(2): ... def K_over_P(T, A=Antoine_As[i], B=Antoine_Bs[i], C=Antoine_Cs[i], fl=fugacities_liq[i], ... fg=fugacities_gas[i], gamma=gammas[i], poy=Poyntings[i]): ... return Antoine(T, A, B, C)*gamma*poy*fl/fg ... funcs.append(K_over_P) >>> _, _, VF, xs, ys = flash_ideal(zs, funcs, Tcs=Tcs, P=1e5, T=364.0) >>> VF, xs, ys (0.510863971792927, [0.5573493403937615, 0.4426506596062385], [0.4450898279593881, 0.5549101720406119]) Note that while this works for PT composition independent flashes - an outer iterating loop is needed for composition dependence! ''' T_MAX = 50000.0 N = len(zs) cmps = range(N) if T is not None and P is not None: P_inv = 1.0 / P Ks = [0.0] * N for i in cmps: Ks[i] = P_inv * funcs[i](T) ans = (T, P) + flash_inner_loop(zs=zs, Ks=Ks) return ans if T is not None and VF == 0.0: ys = [0.0] * N P_bubble = 0.0 for i in cmps: v = funcs[i](T) * zs[i] P_bubble += v ys[i] = v P_inv = 1.0 / P_bubble for i in cmps: ys[i] *= P_inv return (T, P_bubble, 0.0, zs, ys) if T is not None and VF == 1.0: xs = [0.0] * N P_dew = 0. for i in cmps: v = zs[i] / funcs[i](T) P_dew += v xs[i] = v P_dew = 1. / P_dew for i in cmps: xs[i] *= P_dew return (T, P_dew, 1.0, xs, zs) elif T is not None and VF is not None: # Solve for in the middle of Pdew P_low = flash_ideal(zs, funcs, Tcs, T=T, VF=1)[1] P_high = flash_ideal(zs, funcs, Tcs, T=T, VF=0)[1] info = [] def to_solve(P, info): T_calc, P_calc, VF_calc, xs, ys = flash_ideal(zs, funcs, Tcs, T=T, P=P) info[:] = T_calc, P_calc, VF_calc, xs, ys err = VF_calc - VF return err P = brenth(to_solve, P_low, P_high, args=(info, )) return tuple(info) if Tcs is None: # numba: delete Tcs = [fi.solve_prop(1e6) for fi in funcs] # numba: delete if P is not None and VF == 1: def to_solve(T_guess): T_guess = abs(T_guess) P_dew = 0. for i in cmps: P_dew += zs[i] / funcs[i](T_guess) P_dew = 1. / P_dew return P_dew - P # 2/3 average critical point T_guess = .66666 * sum([Tcs[i] * zs[i] for i in cmps]) try: T_dew = abs(secant(to_solve, T_guess, xtol=1e-12, maxiter=50)) except Exception as e: T_dew = None if T_dew is None or T_dew > T_MAX * 5.0: # Went insanely high T, bound it with brenth T_low_guess = sum([.1 * Tcs[i] * zs[i] for i in cmps]) bound = True try: err_low = to_solve(T_low_guess) except: bound = False try: err_high = to_solve(T_MAX) except: bound = False if bound and err_low * err_high > 0.0: bound = False if bound: T_dew = brenth(to_solve, T_low_guess, T_MAX, fa=err_low, fb=err_high) else: T_dew = secant(to_solve, min(min(Tcs) * 0.9, T_guess), xtol=1e-12, maxiter=50, bisection=True, high=min(Tcs)) xs = [P] * N for i in range(N): xs[i] *= zs[i] / funcs[i](T_dew) return (T_dew, P, 1.0, xs, zs) elif P is not None and VF == 0: def to_solve(T_guess): # T_guess = abs(T_guess) P_bubble = 0.0 for i in cmps: P_bubble += zs[i] * funcs[i](T_guess) return P_bubble - P # 2/3 average critical point T_guess = sum([.55 * Tcs[i] * zs[i] for i in cmps]) try: T_bubble = abs( secant(to_solve, T_guess, maxiter=50, bisection=True, xtol=1e-12)) except: T_bubble = None if T_bubble is None or T_bubble > T_MAX * 5.0: # Went insanely high T, bound it with brenth T_low_guess = sum([.1 * Tcs[i] * zs[i] for i in cmps]) bound = True try: err_low = to_solve(T_low_guess) except: bound = False try: err_high = to_solve(T_MAX) except: bound = False if bound and err_low * err_high > 0.0: bound = False if bound: T_bubble = brenth(to_solve, T_low_guess, T_MAX, fa=err_low, fb=err_high) else: Tc_min = min(Tcs) T_bubble = secant(to_solve, min(Tc_min * 0.9, T_guess), maxiter=50, bisection=True, high=Tc_min, xtol=1e-12) P_inv = 1.0 / P ys = [0.0] * N for i in range(N): ys[i] = zs[i] * P_inv * funcs[i](T_bubble) return (T_bubble, P, 0.0, zs, ys) elif P is not None and VF is not None: bound = True try: T_low = flash_ideal(zs, funcs, Tcs, P=P, VF=1)[0] T_high = flash_ideal(zs, funcs, Tcs, P=P, VF=0)[0] except: bound = False info = [] def err(T, zs, funcs, Tcs, P, VF, info, ignore_err): try: T_calc, P_calc, VF_calc, xs, ys = flash_ideal(zs, funcs, Tcs, T=T, P=P) except: if ignore_err: return -0.5 else: raise ValueError("No solution in inner loop") info[:] = T_calc, P_calc, VF_calc, xs, ys return VF_calc - VF if bound: P = brenth(err, T_low, T_high, xtol=1e-14, args=(zs, funcs, Tcs, P, VF, info, False)) else: T_guess = .5 * sum([Tcs[i] * zs[i] for i in cmps]) Tc_min = min(Tcs) # Starting at the lowest component's Tc should guarantee starting at two phases P = secant(err, Tc_min * (1.0 - 1e-7), xtol=1e-12, high=Tc_min, bisection=True, args=(zs, funcs, Tcs, P, VF, info, True)) return tuple(info) else: raise ValueError("Provide two of P, T, and VF")
def flash_Tb_Tc_Pc(zs, Tbs, Tcs, Pcs, T=None, P=None, VF=None): r'''PVT flash model using a model published in [1]_, which provides a PT surface using only each compound's boiling temperature and critical temperature and pressure. This is useful for obtaining initial guesses for more rigorous models, or it can be used as its own model. Capable of solving with two of `T`, `P`, and `VF` for the other one; that results in three solve modes, but for `VF=1` and `VF=0`, there are additional solvers; for a total of seven solvers implemented. This model uses `flash_inner_loop` to solve the Rachford-Rice problem. .. math:: K_i = \frac{P_{c,i}^{\left(\frac{1}{T} - \frac{1}{T_{b,i}} \right) / \left(\frac{1}{T_{c,i}} - \frac{1}{T_{b,i}} \right)}}{P} Parameters ---------- zs : list[float] Mole fractions of the phase being flashed, [-] Tbs : list[float] Boiling temperatures of all species, [K] Tcs : list[float] Critical temperatures of all species, [K] Pcs : list[float] Critical pressures of all species, [Pa] T : float, optional Temperature, [K] P : float, optional Pressure, [Pa] VF : float, optional Molar vapor fraction, [-] Returns ------- T : float Temperature, [K] P : float Pressure, [Pa] VF : float Molar vapor fraction, [-] xs : list[float] Mole fractions of liquid phase, [-] ys : list[float] Mole fractions of vapor phase, [-] Notes ----- For the cases where `VF` is 1 or 0 and T is known, an explicit solution is used. For the same cases where `P` and `VF` are known, there is no explicit solution available. There is an internal `Tmax` parameter, set to 50000 K; which, in the event of convergence of the Secant method, is used as a bounded for a bounded solver. It is used in the PVF solvers. This typically allows pressures up to 2 MPa to be converged to. Failures may still occur for other conditions. This model is based on [1]_, which aims to estimate dew and bubble points using the same K value formulation as used here. While this implementation uses a numerical solver to provide an exact bubble/dew point estimate, [1]_ suggests a sequential substitution and flowchart based solver with loose tolerances. That model was also implemented, but found to be slower and less reliable than this implementation. Examples -------- >>> Tcs = [305.322, 540.13] >>> Pcs = [4872200.0, 2736000.0] >>> Tbs = [184.55, 371.53] >>> zs = [0.4, 0.6] >>> flash_Tb_Tc_Pc(zs=zs, Tcs=Tcs, Pcs=Pcs, Tbs=Tbs, T=300, P=1e5) (300, 100000.0, 0.38070407481453833, [0.03115784303656836, 0.9688421569634316], [0.9999999998827086, 1.172914188751506e-10]) References ---------- .. [1] Kandula, Vamshi Krishna, John C. Telotte, and F. Carl Knopf. "It’s Not as Easy as It Looks: Revisiting Peng—Robinson Equation of State Convergence Issues for Dew Point, Bubble Point and Flash Calculations." International Journal of Mechanical Engineering Education 41, no. 3 (July 1, 2013): 188-202. https://doi.org/10.7227/IJMEE.41.3.2. ''' T_MAX = 50000 N = len(zs) cmps = range(N) # Assume T and P to begin with if T is not None and P is not None: Ks = [ Pcs[i]**((1.0 / T - 1.0 / Tbs[i]) / (1.0 / Tcs[i] - 1.0 / Tbs[i])) / P for i in cmps ] return (T, P) + flash_inner_loop(zs=zs, Ks=Ks, check=True) if T is not None and VF == 0: P_bubble = 0.0 for i in cmps: P_bubble += zs[i] * Pcs[i]**((1.0 / T - 1.0 / Tbs[i]) / (1.0 / Tcs[i] - 1.0 / Tbs[i])) return flash_Tb_Tc_Pc(zs, Tbs, Tcs, Pcs, T=T, P=P_bubble) if T is not None and VF == 1: # Checked to be working vs. PT implementation. P_dew = 0. for i in cmps: P_dew += zs[i] / (Pcs[i]**((1.0 / T - 1.0 / Tbs[i]) / (1.0 / Tcs[i] - 1.0 / Tbs[i]))) P_dew = 1. / P_dew return flash_Tb_Tc_Pc(zs, Tbs, Tcs, Pcs, T=T, P=P_dew) elif T is not None and VF is not None: # Solve for in the middle of Pdew P_low = flash_Tb_Tc_Pc(zs, Tbs, Tcs, Pcs, T=T, VF=1)[1] P_high = flash_Tb_Tc_Pc(zs, Tbs, Tcs, Pcs, T=T, VF=0)[1] info = [] def err(P): T_calc, P_calc, VF_calc, xs, ys = flash_Tb_Tc_Pc(zs, Tbs, Tcs, Pcs, T=T, P=P) info[:] = T_calc, P_calc, VF_calc, xs, ys return VF_calc - VF P = brenth(err, P_low, P_high) return tuple(info) elif P is not None and VF == 1: checker = oscillation_checker() def to_solve(T_guess): T_guess = abs(T_guess) P_dew = 0. for i in range(len(zs)): P_dew += zs[i] / (Pcs[i]**((1.0 / T_guess - 1.0 / Tbs[i]) / (1.0 / Tcs[i] - 1.0 / Tbs[i]))) P_dew = 1. / P_dew err = P_dew - P if checker(T_guess, err): raise ValueError("Oscillation") # print(T_guess, err) return err Tc_pseudo = sum([Tcs[i] * zs[i] for i in cmps]) T_guess = 0.666 * Tc_pseudo try: T_dew = abs(secant(to_solve, T_guess, maxiter=50, ytol=1e-2)) # , high=Tc_pseudo*3 except: T_dew = None if T_dew is None or T_dew > T_MAX * 5.0: # Went insanely high T, bound it with brenth T_low_guess = sum([.1 * Tcs[i] * zs[i] for i in cmps]) checker = oscillation_checker(both_sides=True, minimum_progress=.05) try: T_dew = brenth(to_solve, T_MAX, T_low_guess) except NotBoundedError: raise Exception( "Bisecting solver could not find a solution between %g K and %g K" % (T_MAX, T_low_guess)) return flash_Tb_Tc_Pc(zs, Tbs, Tcs, Pcs, T=T_dew, P=P) elif P is not None and VF == 0: checker = oscillation_checker() def to_solve(T_guess): T_guess = abs(T_guess) P_bubble = 0.0 for i in cmps: P_bubble += zs[i] * Pcs[i]**((1.0 / T_guess - 1.0 / Tbs[i]) / (1.0 / Tcs[i] - 1.0 / Tbs[i])) err = P_bubble - P if checker(T_guess, err): raise ValueError("Oscillation") # print(T_guess, err) return err # 2/3 average critical point Tc_pseudo = sum([Tcs[i] * zs[i] for i in cmps]) T_guess = 0.55 * Tc_pseudo try: T_bubble = abs(secant(to_solve, T_guess, maxiter=50, ytol=1e-2)) # , high=Tc_pseudo*4 except Exception as e: # print(e) checker = oscillation_checker(both_sides=True, minimum_progress=.05) T_bubble = None if T_bubble is None or T_bubble > T_MAX * 5.0: # Went insanely high T (or could not converge because went too high), bound it with brenth T_low_guess = 0.1 * Tc_pseudo try: T_bubble = brenth(to_solve, T_MAX, T_low_guess) except NotBoundedError: raise Exception( "Bisecting solver could not find a solution between %g K and %g K" % (T_MAX, T_low_guess)) return flash_Tb_Tc_Pc(zs, Tbs, Tcs, Pcs, T=T_bubble, P=P) elif P is not None and VF is not None: T_low = flash_Tb_Tc_Pc(zs, Tbs, Tcs, Pcs, P=P, VF=1)[0] T_high = flash_Tb_Tc_Pc(zs, Tbs, Tcs, Pcs, P=P, VF=0)[0] info = [] def err(T): T_calc, P_calc, VF_calc, xs, ys = flash_Tb_Tc_Pc(zs, Tbs, Tcs, Pcs, T=T, P=P) info[:] = T_calc, P_calc, VF_calc, xs, ys return VF_calc - VF P = brenth(err, T_low, T_high) return tuple(info) else: raise ValueError("Provide two of P, T, and VF")