A Python library for simulating one-dimensional heat conduction and exergy conduction through multi-layer solid materials. This library provides tools for analyzing thermal performance and exergy consumption in solid materials using both unsteady-state and pseudo-unsteady-state methods.
The exergy-conduction library enables detailed thermal and exergy analysis of multi-layer solid assemblies. It solves transient heat conduction problems using the Tri-Diagonal Matrix Algorithm (TDMA) and calculates exergy consumption rates throughout the material layers.
-
Unsteady-State Exergy Conduction: Solves transient heat conduction equations considering thermal mass effects, providing time-dependent temperature and exergy consumption profiles.
-
Pseudo-Unsteady-State Exergy Conduction: Assumes steady-state conditions at each time step, calculating exergy consumption based on instantaneous temperature gradients without considering thermal mass.
- Multi-Layer Material Modeling: Define complex multi-layer assemblies with multiple material layers
- Unsteady-State Simulation: Transient heat conduction analysis with thermal mass effects
- Pseudo-Unsteady-State Simulation: Steady-state analysis at each time step
- Flexible Boundary Conditions: Support for fixed temperature and convection boundary conditions
- TDMA Solver: Efficient tri-diagonal matrix algorithm for solving transient heat conduction
- Comprehensive Results: Temperature profiles, heat fluxes, and exergy consumption rates
- Built-in Unit Conversions: Temperature conversion utilities and unit constants
- Python >= 3.11
- NumPy >= 2.3.2
- Pandas >= 2.3.1
- tqdm >= 4.67.1
pip install numpy pandas tqdmOr install from the project directory:
pip install -e .import numpy as np
from exergy_conduction.main import Layer, Construction
# Define material properties for concrete
concrete = Layer(
L=0.2, # Thickness: 20 cm
dx=0.01, # Discretization: 1 cm
k=1.4, # Thermal conductivity [W/mK]
c=880, # Specific heat capacity [J/kgK]
rho=2400 # Density [kg/m³]
)
# Create multi-layer assembly
assembly = Construction(layers=[concrete])
# Define boundary conditions (24 hours, hourly data)
hours = 24
T_left = np.array([5.0] * hours) # Left boundary temperature [°C]
T_right = np.array([20.0] * hours) # Right boundary temperature [°C]
T_env = np.array([10.0] * hours) # Environment temperature for exergy [°C]
T_initial = 12.5 # Initial temperature [°C]
dt = 3600 # Time step: 1 hour [s]
# Run unsteady-state simulation
assembly.solve_unsteady_exergy_conduction(
T0=T_env,
Left_BC=T_left,
Right_BC=T_right,
T_init=T_initial,
dt=dt
)
# Access results
results = assembly.results['unsteady']
temperatures = results['T'] # Node temperatures [K]
heat_flux = results['q'] # Heat flux [W/m²]
exergy_consumption = results['CXcR'] # Exergy consumption rate [W/m²]
# Print simulation summary
assembly.print_simulation_summary('unsteady')import numpy as np
from exergy_conduction.main import Layer, Construction
# Define multiple layers
layer1 = Layer(L=0.02, dx=0.01, k=0.7, c=1000, rho=1800)
layer2 = Layer(L=0.10, dx=0.01, k=0.04, c=1030, rho=30)
layer3 = Layer(L=0.15, dx=0.01, k=1.4, c=880, rho=2400)
layer4 = Layer(L=0.02, dx=0.01, k=0.7, c=1000, rho=1800)
# Create multi-layer assembly
assembly = Construction(layers=[
layer1,
layer2,
layer3,
layer4
])
# Set convection boundary conditions
assembly.boundary_layer(
h_c_left=25.0, # Left boundary convection coefficient [W/m²K]
h_c_right=8.0 # Right boundary convection coefficient [W/m²K]
)
# Define time-varying boundary conditions
hours = 24
time = np.arange(hours)
T_left = 10 + 5 * np.sin(2 * np.pi * time / 24) # Sinusoidal left boundary temp [°C]
T_right = np.array([22.0] * hours) # Constant right boundary temp [°C]
T_env = np.array([12.0] * hours) # Environment temp [°C]
T_initial = 15.0 # Initial temp [°C]
dt = 3600 # 1 hour time step [s]
# Run simulation
assembly.solve_unsteady_exergy_conduction(
T0=T_env,
Left_BC=T_left,
Right_BC=T_right,
T_init=T_initial,
dt=dt
)
# Access layer-specific results
print(f"Total thickness: {assembly.thick:.3f} m")
print(f"Total thermal resistance: {assembly.R_tot:.4f} m²K/W")
print(f"Number of discretization nodes: {assembly.DN}")import numpy as np
from exergy_conduction.main import Layer, Construction
# Define material layer
concrete = Layer(L=0.2, dx=0.01, k=1.4, c=880, rho=2400)
assembly = Construction(layers=[concrete])
# Set convection boundaries
assembly.boundary_layer(h_c_left=25.0, h_c_right=8.0)
# Define boundary conditions
hours = 24
T_left = 10 + 5 * np.sin(2 * np.pi * np.arange(hours) / 24)
T_right = np.array([22.0] * hours)
T_env = np.array([12.0] * hours)
# Run pseudo-unsteady simulation
assembly.solve_pseudo_unsteady_exergy_conduction(
T0=T_env,
Left_BC=T_left,
Right_BC=T_right
)
# Access results
results = assembly.results['pseudo_unsteady']
temperatures = results['T'] # Node temperatures [K]
interface_temps = results['T_intf_BC'] # Interface temperatures [K]
heat_flux = results['q'] # Heat flux [W/m²]
exergy_consumption = results['CXcR'] # Exergy consumption rate [W/m²]
# Compare with unsteady results
print("Pseudo-unsteady simulation completed")
print(f"Average heat flux: {heat_flux.mean():.2f} W/m²")
print(f"Total exergy consumption: {exergy_consumption.sum():.6f} W/m²")Convert temperature from Celsius to Kelvin.
Parameters:
temp_C: Temperature in Celsius [°C]
Returns:
- Temperature in Kelvin [K]
Convert temperature from Kelvin to Celsius.
Parameters:
temp_K: Temperature in Kelvin [K]
Returns:
- Temperature in Celsius [°C]
Calculate midpoint values between consecutive elements along specified axis.
Parameters:
mat: Input array (1D or 2D)axis: Axis along which to interpolate (0: row direction, 1: column direction)
Returns:
- Array with midpoint values (one element smaller along specified axis)
Remove the last row from an array. Used for adjusting results when extra time steps are calculated.
Represents a single material layer in a multi-layer solid assembly.
L(float): Layer thickness [m]dx(float): Spatial discretization size [m]k(float): Thermal conductivity [W/mK]c(float): Specific heat capacity [J/kgK]rho(float): Density [kg/m³]
DN(int): Number of discretization nodesC(float): Volumetric heat capacity [J/m³K] =c * rhoR(float): Thermal resistance [m²K/W] =dx / kK(float): Thermal conductance [W/m²K] =k / dxalpha(float): Thermal diffusivity [m²/s] =k / C
layer = Layer(L=0.2, dx=0.01, k=1.4, c=880, rho=2400)
print(f"Thermal resistance: {layer.R:.4f} m²K/W")
print(f"Number of nodes: {layer.DN}")Represents a multi-layer solid assembly with methods for thermal and exergy analysis.
Construction(layers: List[Layer])Parameters:
layers: List ofLayerobjects in order from left to right
layer_num(int): Number of layersDN(int): Total number of discretization nodesthick(float): Total thickness [m]R_tot(float): Total thermal resistance [m²K/W]K_tot(float): Total thermal conductance [W/m²K]node_pos(np.ndarray): Node center positions from left surface [cm]node_pos_with_BC(np.ndarray): Node positions including boundaries [cm]cell_surf_pos(np.ndarray): Cell surface positions [cm]
Set convection heat transfer coefficients at boundaries.
Parameters:
h_c_left(float, optional): Left boundary convection coefficient [W/m²K]. Default:None(infinite, fixed temperature)h_c_right(float, optional): Right boundary convection coefficient [W/m²K]. Default:None(infinite, fixed temperature)
Note: When convection coefficients are set, boundary conditions represent fluid/air node temperatures rather than surface temperatures. This is applicable to any convective heat transfer scenario, not limited to air.
Solve unsteady-state exergy conduction problem.
Parameters:
T0(np.ndarray): Environment temperature array [°C], shape(tN,)Left_BC(np.ndarray): Left boundary temperature array [°C], shape(tN,)Right_BC(np.ndarray): Right boundary temperature array [°C], shape(tN,)T_init(float): Initial temperature [°C]dt(float): Time step [s]
Returns:
- None (results stored in
self.results['unsteady'])
Results Dictionary:
T: Node temperatures [K], shape(tN-1, DN)T_L: Left surface temperatures [K], shape(tN-1, DN)T_R: Right surface temperatures [K], shape(tN-1, DN)q: Average heat flux [W/m²], shape(tN-1, DN)q_in: Incoming heat flux [W/m²], shape(tN-1, DN)q_out: Outgoing heat flux [W/m²], shape(tN-1, DN)T_half: Half-time node temperatures [K], shape(tN-1, DN)T_L_half: Half-time left surface temperatures [K], shape(tN-1, DN)T_R_half: Half-time right surface temperatures [K], shape(tN-1, DN)T0_half: Half-time environment temperatures [K], shape(tN-1, 1)q_half: Half-time heat flux [W/m²], shape(tN-1, DN)CXcR: Exergy consumption rate [W/m²], shape(tN-1, DN)dt: Time step [s]simulation_type: 'unsteady'
Note: Results arrays have tN-1 rows because the last time step is eliminated (see eliminate_last_row function).
Solve pseudo-unsteady-state exergy conduction problem.
Parameters:
T0(np.ndarray): Environment temperature array [°C], shape(tN,)Left_BC(np.ndarray): Left boundary temperature array [°C], shape(tN,)Right_BC(np.ndarray): Right boundary temperature array [°C], shape(tN,)
Returns:
- None (results stored in
self.results['pseudo_unsteady'])
Results Dictionary:
T: Node temperatures [K], shape(tN, DN)T_intf: Interface temperatures [K], shape(tN, DN)T_surf_left: Left surface temperatures [K], shape(tN, 1)T_surf_right: Right surface temperatures [K], shape(tN, 1)T_intf_BC: Interface temperatures with boundaries [K], shape(tN, DN+1)q: Heat flux [W/m²], shape(tN, 1)(spatially constant)q_intf: Interface heat flux [W/m²], shape(tN, DN+1)q_node: Node heat flux [W/m²], shape(tN, DN)CXcR: Exergy consumption rate [W/m²], shape(tN, DN)R_L_cumsum: Cumulative left resistancesR_cumsum: Cumulative resistancessimulation_type: 'pseudo_unsteady'
Print comprehensive simulation summary.
Parameters:
simulation_type(str): Type of simulation ('unsteady' or 'pseudo_unsteady')
Output includes:
- Construction information (layers, thickness, nodes, thermal resistance)
- Simulation information (time steps, duration)
- Boundary condition information
- Temperature analysis (ranges, variations)
- Heat flux analysis
- Exergy consumption analysis
- Layer details
Get results dictionary for specified simulation type.
Parameters:
simulation_type(str): Type of simulation ('unsteady' or 'pseudo_unsteady')
Returns:
- Dictionary containing simulation results, or
Noneif not available
List all available simulation results.
Returns:
- None (prints available results to console)
- Input temperatures: Celsius [°C]
- Internal calculations: Kelvin [K]
- Output temperatures: Kelvin [K] (convert using
K2C()if needed)
- Thermal conductivity (
k): [W/mK] - Specific heat capacity (
c): [J/kgK] - Density (
rho): [kg/m³] - Volumetric heat capacity (
C): [J/m³K] - Thermal resistance (
R): [m²K/W] - Thermal conductance (
K): [W/m²K] - Thermal diffusivity (
alpha): [m²/s]
- Thickness (
L): [m] - Discretization (
dx): [m] - Node positions: [cm] (converted from meters)
- Time step (
dt): [s] - Time arrays: Number of time steps
- Heat flux (
q): [W/m²] - Convection coefficient (
h_c): [W/m²K] - Exergy consumption rate (
CXcR): [W/m²]
-
Fixed Temperature Boundaries (default):
- When
boundary_layer()is not called or coefficients areNone, boundaries represent surface temperatures directly. h_c_leftandh_c_rightare set to infinity.
- When
-
Convection Boundaries:
- When convection coefficients are specified, boundary conditions represent fluid/air node temperatures.
- Surface temperatures are calculated considering convection resistance.
- The relationship:
T_surface = f(T_fluid, T_node, R_convection, R_material)
-
Unsteady-State: Considers thermal mass effects. Temperature changes propagate through the material over time. Requires initial temperature and time step. More computationally intensive but physically accurate for transient analysis.
-
Pseudo-Unsteady-State: Assumes steady-state at each time step. No thermal mass effects. Faster computation but less accurate for rapid temperature changes. Suitable for quasi-steady conditions.
- Layer thickness
Lmust be divisible by discretization sizedx. - The library automatically calculates the number of nodes:
DN = L / dx - Each layer can have different
dxvalues.
Results are stored in two ways:
self.results[simulation_type]: Dictionary with all results- Direct attributes:
self.T,self.q,self.CXcR, etc. (for backward compatibility)
The TDMA solver implementation follows:
- Reference: https://doi.org/10.1016/j.ijheatmasstransfer.2017.09.057 [Appendix B - Eq.(B7)]