NonLumpedSlab¶
- class pychemengg.heattransfer.transient.NonLumpedSlab(thickness=None, surfacearea=None, volume=None, density=None, specificheat=None, thermalconductivity=None, thermaldiffusivity=None, heattransfercoefficient=None, T_infinity=None, T_initial=None)[source]¶
Bases:
object
Model for nonlumped analysis of rectangular solid object.
- Parameters
- thicknessint or float
Thickness of solid object
- surfaceareaint or float
Surface area of solid object.
- volumeint or float
Volume of solid object.
- densityint or float
Density of solid object.
- specificheatint or float
Specific heat of solid object.
- thermalconductivityint or float
Thermal conductivity of solid object.
- thermaldiffusivityint or float
Thermal diffusivity of solid object.
- heattransfercoefficientint or float
Heat transfer coefficient between solid object and surrounding.
- T_infinityint or float
Temperature of surroundings.
- T_initialint or float
Temperature of solid object at time = 0.
Examples
First import the module transient
Units used in this example: SI system
However, any consistent units can be used
>>> from pychemengg.heattransfer import transient >>> plate = transient.NonLumpedSlab(thickness=4e-2, surfacearea=1, volume=1*4e-2, density=8530, specificheat=380, thermalconductivity=110, thermaldiffusivity=None, heattransfercoefficient=120, T_infinity=500, T_initial=20) # This will create an instance of 'NonLumpedSlab' with a name 'plate'
- Attributes
- See “Parameters”. All parameters are attributes. Additional attributes are listed below.
- massint or float
Mass of solid object computed as (volume * density) of solid object.
- __init__(thickness=None, surfacearea=None, volume=None, density=None, specificheat=None, thermalconductivity=None, thermaldiffusivity=None, heattransfercoefficient=None, T_infinity=None, T_initial=None)[source]¶
Methods
__init__
([thickness, surfacearea, volume, …])calc_Bi
()Computes Biot number.
calc_Fo
([time])Computes Fourier number.
calc_eigenvalues
([numberof_eigenvalues_desired])Computes eigen values of characteristic equation for Slab geometry.
calc_heatrateof_conv_at_time_t
([time])Heat rate of convection between object and surroundings at a given time = t.
Maximum possible heat transfer between solid object and surroundings.
calc_temperature_of_solid_at_time_t
([time, …])Calculates temperature of solid object at a given time = t and position = x.
Heat transferred between solid object and surroundings during time interval = 0 to t.
- calc_Bi()[source]¶
Computes Biot number.
- Parameters
- None_requiredNone
Attributes that are already defined are used in calculation.
- Returns
- Biint or float
Biot number
Notes
Biot number is calculated using the following formula.
\[Bi = \frac {h L_{c}} {k} \]where:
h = heat transfer coefficient
k = thermal conductivity of solid object
\(L_c\) = characteristic length of solid object = thickness/2
Bi = Biot number
References
[1] Y. A. Cengel and A. J. Ghajar, “Heat And Mass Transfer Fundamentals and Applications”, 6th Edition. New York, McGraw Hill Education, 2020.
Examples
First import the module transient
Units used in this example: SI system
However, any consistent units can be used
>>> from pychemengg.heattransfer import transient >>> plate = transient.NonLumpedSlab(thickness=4e-2, surfacearea=1, volume=1*4e-2, density=8530, specificheat=380, thermalconductivity=110, thermaldiffusivity=None, heattransfercoefficient=120, T_infinity=500, T_initial=20) # This will create an instance of 'NonLumpedSlab' with a name 'plate' # Next call calc_Bi >>> plate.calc_Bi() 0.021818181818181816
- calc_Fo(time=None)[source]¶
Computes Fourier number.
- Parameters
- timeint or float
Time at which temperature or heat transfer is to be evaluated.
- Returns
- Foint or float
Fourier number
Notes
Fourier number is calculated using the following formula.
\[Fo = \frac {\alpha t} {L_c^2} \]where:
\(\alpha\) = thermal diffusivity
t = time at which temperature or heat transfer is to be evaluated
\(L_c\) = characteristic length = (slab thickness)/2
Fo = Fourier number
References
[1] Y. A. Cengel and A. J. Ghajar, “Heat And Mass Transfer Fundamentals and Applications”, 6th Edition. New York, McGraw Hill Education, 2020.
Examples
First import the module transient
Units used in this example: SI system
However, any consistent units can be used
>>> from pychemengg.heattransfer import transient >>> plate = transient.NonLumpedSlab(thickness=4e-2, surfacearea=1, volume=1*4e-2, density=8530, specificheat=380, thermalconductivity=110, thermaldiffusivity=None, heattransfercoefficient=120, T_infinity=500, T_initial=20) # This will create an instance of 'NonLumpedSlab' with a name 'plate' # Next call calc_Fo assuming temperature is required at 7 min >>> plate.calc_Fo(time=7*60) 35.63275128031097
- calc_eigenvalues(numberof_eigenvalues_desired=10)[source]¶
Computes eigen values of characteristic equation for Slab geometry.
- Parameters
- numberof_eigenvalues_desiredint or float (default = 10)
Number of eigen values desired for the characteristic equation.
- Returns
- eigenvaluesnp.array of int or float
Eigen values
Notes
Eigen values are calculated as roots of the following equation.
\[x_n tan(x_n) - Bi = 0 , n = 1 \hspace{2pt} to \hspace{2pt} \infty \]where:
\(x_n\) = nth eigen value
Bi = Biot number
References
[1] Y. A. Cengel and A. J. Ghajar, “Heat And Mass Transfer Fundamentals and Applications”, 6th Edition. New York, McGraw Hill Education, 2020.
Examples
First import the module transient
Units used in this example: SI system
However, any consistent units can be used
>>> from pychemengg.heattransfer import transient >>> plate = transient.NonLumpedSlab(thickness=4e-2, surfacearea=1, volume=1*4e-2, density=8530, specificheat=380, thermalconductivity=110, thermaldiffusivity=None, heattransfercoefficient=120, T_infinity=500, T_initial=20) # This will create an instance of 'NonLumpedSlab' with a name 'plate' # Next call calc_Bi >>> plate.calc_Bi() 0.021818181818181816 # Let first 5 eigen values be required >>> plate.calc_eigenvalues(numberof_eigenvalues_desired=5) array([ 0.14717481, 3.1485222 , 6.28665585, 9.42709237, 12.56810661])
- calc_heatrateof_conv_at_time_t(time=None)[source]¶
Heat rate of convection between object and surroundings at a given time = t.
- Parameters
- timeint or float
Time instant from begining of process, at which heat rate is to be found.
- Returns
- heat rate of convectionint or float ; Positive: Heat is gained by object, Negative: Heat is lost by object
Heat rate of convection between solid object and surroundings at time = t.
Notes
Heat rate of convection is calculated using the following formula:
\[q_{t} = h A_s (T_{infinity} - T_{t}) \]where:
t = time at which temperature is to be computed
h = heat transfer coefficient
\(T_{infinity}\) = temperature of surrounding fluid
\(T_{t}\) = temperature of surface of solid object at time = t
\(A_s\) = surface area of solid object
\(q_{t}\) = heat rate of convection at time = t
References
[1] Y. A. Cengel and A. J. Ghajar, “Heat And Mass Transfer Fundamentals and Applications”, 6th Edition. New York, McGraw Hill Education, 2020.
Examples
First import the module transient
Units used in this example: SI system
However, any consistent units can be used
>>> from pychemengg.heattransfer import transient >>> plate = transient.NonLumpedSlab(thickness=4e-2, surfacearea=1, volume=1*4e-2, density=8530, specificheat=380, thermalconductivity=110, thermaldiffusivity=None, heattransfercoefficient=120, T_infinity=500, T_initial=20) # This will create an instance of 'NonLumpedSlab' with a name 'plate' # Let temperature at time = 7 min after start of the process be required. # Next call the following >>> plate.calc_Bi() 0.021818181818181816 >>> plate.calc_Fo(time=7*60) 35.63275128031097 plate.calc_eigenvalues() array([ 0.14717481, 3.1485222 , 6.28665585, 9.42709237, 12.56810661, 15.70935213, 18.85071334, 21.99214066, 25.13360932, 28.27510552]) >>> plate.calc_temperature_of_solid_at_time_t(time=7*60, xposition_tofindtemp=plate.thickness/2) 279.76430920417204 # Next call the following >>> plate.calc_heatrateof_conv_at_time_t(time=7*60) 26428.282895499357 # Positive sign indicates gain of heat by the solid object
- calc_maxheattransferpossible()[source]¶
Maximum possible heat transfer between solid object and surroundings.
- Parameters
- None_requiredNone
Attributes that are already defined are used in calculation.
- Returns
- maximum heat transfer possible: int or float; Positive: Heat is gained by object, Negative: Heat is lost by object
Maximum heat transfer posssible between object and surroundings.
Notes
Maximum heat transfer possible between solid object and surroundings is calculated using the following formula. This is based on the assumption that final object temperature will eventually reach surrounding temperature of \(T_{infinity}\)
\[q_{max} = m C_p (T_{infinity} - T_{initial}) \]where:
m = mass of solid object
\(C_{p}\) = specific heat of solid object
\(T_{infinity}\) = temperature of surrounding, which the solid object will eventually attain
\(T_{initial}\) = temperature of solid object at time = initial
\(q_{max}\) = max heat transfer possible
References
[1] Y. A. Cengel and A. J. Ghajar, “Heat And Mass Transfer Fundamentals and Applications”, 6th Edition. New York, McGraw Hill Education, 2020.
Examples
First import the module transient
Units used in this example: SI system
However, any consistent units can be used
>>> from pychemengg.heattransfer import transient >>> plate = transient.NonLumpedSlab(thickness=4e-2, surfacearea=1, volume=1*4e-2, density=8530, specificheat=380, thermalconductivity=110, thermaldiffusivity=None, heattransfercoefficient=120, T_infinity=500, T_initial=20) # This will create an instance of 'NonLumpedSlab' with a name 'plate' >>> plate.calc_maxheattransferpossible() 62234880.0 # negative value indicates heat is being lost by the solid object
- calc_temperature_of_solid_at_time_t(time=None, xposition_tofindtemp=None)[source]¶
Calculates temperature of solid object at a given time = t and position = x.
- Parameters
- timeint or float
Time instant from begining of process, at which temperature of solid object is to be found.
- xposition_tofindtempint or float
Distance measured from center of rectangular object where temperature is to be found.
- Returns
- temperatureint or float
Temperature of solid object at time = t and position = x.
Notes
Temperature of solid object at time = t and position = x is calculated using the following formula:
\[T(t) = T_{infinity} + (T_{initial} - T_{infinity}) \displaystyle\sum_{n=1}^\infty \cfrac{4sin(\lambda_n)}{2 \lambda_n + sin(2 \lambda_n)} e^{- \lambda_n^2 \tau} cos(\lambda_n x/L) \]where:
\(T_{infinity}\) = temperature of surrounding fluid
\(T_{initial}\) = intitial temperature of solid object
\(\lambda_n\) = nth eigen value of \(x_n tan(x_n) - Bi = 0\) , n = 1 \(\hspace{2pt} to \hspace{2pt} \infty\)
\(\tau\) = Fourier number
x = distance from center of solid slab where temperature is required (x = 0 for center of slab)
L = thickness/2
References
[1] Y. A. Cengel and A. J. Ghajar, “Heat And Mass Transfer Fundamentals and Applications”, 6th Edition. New York, McGraw Hill Education, 2020.
Examples
First import the module transient
Units used in this example: SI system
However, any consistent units can be used
>>> from pychemengg.heattransfer import transient >>> plate = transient.NonLumpedSlab(thickness=4e-2, surfacearea=1, volume=1*4e-2, density=8530, specificheat=380, thermalconductivity=110, thermaldiffusivity=None, heattransfercoefficient=120, T_infinity=500, T_initial=20) # This will create an instance of 'NonLumpedSlab' with a name 'plate' # Let temperature at time = 7 min after start of the process be required. # Next call the following >>> plate.calc_Bi() 0.021818181818181816 >>> plate.calc_Fo(time=7*60) 35.63275128031097 plate.calc_eigenvalues() array([ 0.14717481, 3.1485222 , 6.28665585, 9.42709237, 12.56810661, 15.70935213, 18.85071334, 21.99214066, 25.13360932, 28.27510552]) >>> plate.calc_temperature_of_solid_at_time_t(time=7*60, xposition_tofindtemp=plate.thickness/2) 279.76430920417204
- calc_totalheat_transferred_during_interval_t()[source]¶
Heat transferred between solid object and surroundings during time interval = 0 to t.
- Parameters
- None_requiredNone
Attributes that are already defined or calculated are used in calculation.
- Returns
- total heat transferredint or float; Positive: Heat is gained by object, Negative: Heat is lost by object
Total heat transferred between object and surroundings during interval 0 to t
Notes
Total heat transferred in interval 0 to t is calculated using the following formula:
\[q_{0 \to t} = q_{max} \left( 1 - \displaystyle\sum_{n=1}^\infty \cfrac{4sin( \lambda_n)}{2 \lambda_n + sin(2 \lambda_n)} \frac{sin( \lambda_n)}{\lambda_n} e^{- \lambda_n^2 \tau} \right) \]where:
\(\lambda_n\) = nth eigen value of \(x_n tan(x_n) - Bi = 0\) , n = 1 \(\hspace{2pt} to \hspace{2pt} \infty\)
\(\tau\) = Fourier number
\(q_{max}\) = maximum heat transfer possible between object and surroundings
References
[1] Y. A. Cengel and A. J. Ghajar, “Heat And Mass Transfer Fundamentals and Applications”, 6th Edition. New York, McGraw Hill Education, 2020.
Examples
First import the module transient
Units used in this example: SI system
However, any consistent units can be used
>>> from pychemengg.heattransfer import transient >>> plate = transient.NonLumpedSlab(thickness=4e-2, surfacearea=1, volume=1*4e-2, density=8530, specificheat=380, thermalconductivity=110, thermaldiffusivity=None, heattransfercoefficient=120, T_infinity=500, T_initial=20) # This will create an instance of 'NonLumpedSlab' with a name 'plate' # Let temperature at time = 7 min after start of the process be required. # Next call the following >>> plate.calc_Bi() 0.021818181818181816 >>> plate.calc_Fo(time=7*60) 35.63275128031097 plate.calc_eigenvalues() array([ 0.14717481, 3.1485222 , 6.28665585, 9.42709237, 12.56810661, 15.70935213, 18.85071334, 21.99214066, 25.13360932, 28.27510552]) >>> plate.calc_totalheat_transferred_during_interval_t() 33472028.92491645 # Positive value means heat is gained by the object.