bloc.reactors#

Extension of Cantera Reactor classes to include new features.

Implementation is based on the examples given in https://cantera.org/3.1/examples/python/reactors/custom2.html

Attributes#

Classes#

CarbonBlackIdealGasReactor

IdealGasReactor with additional features for carbon quality calculations.

CarbonBlackIdealGasConstPressureReactor

IdealGasConstPressureReactor with additional features for carbon quality calculations.

CarbonBlackIdealGasMoleReactor

IdealGasMoleReactor with additional features for carbon quality calculations.

Functions#

get_H2_yield(states)

Compute the H2 yield of the process from the states object.

find_temperature_from_enthalpy(h_mass, X, P_bar, mechanism)

Return the temperature in °C for a given enthalpy in J/kg and for a given gas composition and pressure.

switch_mechanism(gas, new_mechanism[, htol, Xtol, verbose])

Switch the mechanism of a gas object.

solve_isothermal_reactor(qm_kgs, T_reactor_C, t_res, ...)

Solve the kinetics in an isothermal isobaric reactor at T_reactor_C (°C) and P_bar (bar) for a residence time t_res (s).

solve_adiabatic_reactor(qm_kgs, P_in_kW, t_res, X0, ...)

Solve the kinetics in an adiabatic isobaric reactor at P_bar (bar) for a residence time t_res (s).

mix_two_streams(gas_1, qm_1, gas_2, qm_2)

Mix two streams of gases and return the resulting gas object.

compute_reactor_length(qm_in, S_in, states)

qm_in = v * S_in * density <=> v = qm_in / (S_in*density).

compute_reactor_dimensions_with_form_factor(qm_in, ...)

get_reactor_mass(L_reactor, d_reactor, e_insul_layers, ...)

Compute the mass of the reactor based on its dimensions and the insulation layers.

thermal_resistance_radial_conduction(diameter, length, ...)

Compute the thermal resistance of the radial heat transfer in cylindrical geometry.

thermal_resistance_CC(h_CC, S)

Compute the thermal resistance of the conducto-convective heat transfer.

compute_hCC_natConv_vertCyl(diameter, length, ...[, ...])

Compute the h coefficient for natural convection at the external wall of a vertical cylinder.

compute_hrad_linearRadiation(eps, Tmoy)

compute_heat_losses_linear(Tr_C, Tamb_C, R_list)

Solve the heat transfer problem in the linear approximation.

compute_gas2wall_radiative_flux(T_g, T_w, kappa_grey, D)

Compute the radiative flux from the gas to the wall in W/m2.

compute_recovered_power(qm_torch, T_torch_input_C, ...)

Compute the power recovered from the preheating of the torch and the second injection.

compute_residual_heat(P_recovered, gas, qm_tot, T_amb_C)

Module Contents#

class bloc.reactors.CarbonBlackIdealGasReactor(*args, **kwargs)#

Bases: cantera.ExtensibleIdealGasReactor

IdealGasReactor with additional features for carbon quality calculations.

This class extends the cantera.IdealGasReactor class to include additional features for carbon quality calculations. The class is intended to be used in Reactor Networks.

Examples

residence_time#

Residence time in (s). Residence time must be set by the user during simulation.

after_initialize(t0)#
guess_specific_surface(Model=KirkOthmer2004SurfaceArea)#

Evaluate the specific surface of the carbon based on temperature and residence time.

Note

The residence time must be calculated before calling this method; for instance with:

reactor.residence_time = reactor.volume / flow_rate

Or:

for t in logspace(0.001, 1, 50):  # s
    sim.advance(t)
    reactor.residence_time = t
Parameters:

Model (class) – Surface area model to use. Default is KirkOthmer2004SurfaceArea.

Examples

class bloc.reactors.CarbonBlackIdealGasConstPressureReactor(*args, **kwargs)#

Bases: cantera.ExtensibleIdealGasConstPressureReactor

IdealGasConstPressureReactor with additional features for carbon quality calculations.

This class extends the cantera.IdealGasReactor class to include additional features for carbon quality calculations. The class is intended to be used in Reactor Networks.

Examples

Carbon Black Reactor.

Carbon Black Reactor.
residence_time#

Residence time in (s). Residence time must be set by the user during simulation.

after_initialize(t0)#
guess_specific_surface(Model=KirkOthmer2004SurfaceArea)#

Evaluate the specific surface of the carbon based on temperature and residence time.

Note

The residence time must be calculated before calling this method; for instance with:

reactor.residence_time = reactor.volume / flow_rate

Or:

for t in logspace(0.001, 1, 50):  # s
    sim.advance(t)
    reactor.residence_time = t
Parameters:

Model (class) – Surface area model to use. Default is KirkOthmer2004SurfaceArea.

Examples

Carbon Black Reactor.

Carbon Black Reactor.
class bloc.reactors.CarbonBlackIdealGasMoleReactor(*args, **kwargs)#

Bases: cantera.ExtensibleIdealGasMoleReactor

IdealGasMoleReactor with additional features for carbon quality calculations.

This class extends the cantera.IdealGasMoleReactor class to include additional features for carbon quality calculations. The class is intended to be used in Reactor Networks.

Examples

residence_time#

Residence time in (s). Residence time must be set by the user during simulation.

after_initialize(t0)#
guess_specific_surface(Model=KirkOthmer2004SurfaceArea)#

Evaluate the specific surface of the carbon based on temperature and residence time.

Note

The residence time must be calculated before calling this method; for instance with:

reactor.residence_time = reactor.volume / flow_rate

Or:

for t in logspace(0.001, 1, 50):  # s
    sim.advance(t)
    reactor.residence_time = t
Parameters:

Model (class) – Surface area model to use. Default is KirkOthmer2004SurfaceArea.

Examples

bloc.reactors.get_H2_yield(states)#

Compute the H2 yield of the process from the states object.

The H2 yield is define as: amount of H2 produced / maximum amount of H2 that could have been produced. The initial amount of H2 in the feedstock must be subtracted, because it is not produced. Thus: H2_yield = (Y_H2_final - Y_H2_initial) / (Y_H_tot - Y_prod_i)

bloc.reactors.find_temperature_from_enthalpy(h_mass, X, P_bar, mechanism)#

Return the temperature in °C for a given enthalpy in J/kg and for a given gas composition and pressure.

Parameters:
  • h_mass (-) – Specific enthalpy target in J/kg

  • X (-) – Composition of the gas in mole fraction

  • P_bar (-) – Pressure in bar

  • mechanism (-) – Path to the mechanism file

bloc.reactors.switch_mechanism(gas, new_mechanism, htol=0.0001, Xtol=0.0001, verbose=False)#

Switch the mechanism of a gas object.

Useful to switch from plasma to reactor simulation for instance. As some species may not be present in the new mechanism, we compute the adjust the temperature of the new gas object to match the enthalpy of the old gas object, so that energy is conserved. Still, if the chemical enthalpy variation exceed htol, an error is raised. Similarly, if the mole fraction variation exceed Xtol, a error is raised.

Parameters:
  • gas (ct.Solution) – Gas object

  • new_mechanism (str) – Path to the new mechanism file

  • htol (float) – Tolerance for the chemical enthalpy variation between the two mechanisms. Default is 1e-4.

  • Xtol (float) – Tolerance for the mole fraction variation between the two mechanisms. Default is 1e-4.

  • verbose (bool) – If True, print information about the process. Default is False.

Returns:

gas_new – Gas object with the new mechanism

Return type:

ct.Solution

bloc.reactors.solve_isothermal_reactor(qm_kgs, T_reactor_C, t_res, X0, T0_C, P_bar, mechanism=None, dt_short=1e-09, rtolT_prevention=0.0001, rtolT_fatal=0.01, heat_recycling=True, conversion_target=None, H2_yield_target=None, name=None, verbose=2, termination='residence_time', L_r=None, S_r=None)#

Solve the kinetics in an isothermal isobaric reactor at T_reactor_C (°C) and P_bar (bar) for a residence time t_res (s).

Parameters:
  • qm_kgs (-) – Mass flow rate of the input gas in kg/s

  • T_reactor_C (-) – Temperature of the reactor in °C

  • t_res (-) – Residence time in the reactor in s

  • X0 (-) – Composition of the input gas in mole fraction

  • T0_C (-) – Temperature of the input gas in °C

  • P_bar (-) – Pressure in the reactor in bar

  • mechanism (-) – Path to the mechanism file

  • dt_short (-) – Short time step for the simulation. A short time step is used to prevent temperature variations.

  • rtolT_prevention (-) – Relative tolerance for temperature variations. If the relative variation of the temperature exceed rtolT_prevention, the time step is reduced.

  • rtolT_fatal (-) – Relative tolerance for temperature variations. If the relative variation of the temperature exceed rtolT_fatal, the simulation is stopped.

  • heat_recycling (-) – If True, the energy released in exothermic reactions is used in later endothermic reactions. Default is True.

  • conversion_target (-) – Target for the conversion of the input gas. Default is None.

  • H2_yield_target (-) – Target for the H2 yield. Default is None.

  • name (-) – Name of the reactor. Default is None.

  • verbose (int) – if >0, print infos. If >=2, add calculation progress bar.

Notes

Energy calculation: Two terms are computed: - Ecost_init: the energy required to heat the input gas from T0_C to T_reactor_C. Multiplying it by the mass flow

yields the heating power

  • delta_h: the energy required to compensate for the endothermic reactions and maintain the temperature constant in the reactor.

    If heat_recycling is True, the energy released in exothermic reactions is used in later endothermic reactions. Multiplying it by the mass flow yields the chemical power.

To compute delta_h, we use ‘energy = on’ in the reactor object definition, and evaluate the enthalpy variation at each time step. The time step is choosed so that the relative variation of the temperature is less than rtolT_prevention.

bloc.reactors.solve_adiabatic_reactor(qm_kgs, P_in_kW, t_res, X0, T0_C, P_bar, mechanism, dt_short=1e-09, conversion_target=None, H2_yield_target=None, termination='residence_time', L_reactor=0.0, S_reactor=1.0, name=None, verbose=2)#

Solve the kinetics in an adiabatic isobaric reactor at P_bar (bar) for a residence time t_res (s).

Parameters:
  • qm_kgs (-) – Mass flow rate of the input gas in kg/s

  • P_in_kW (-) – Power input in kW

  • t_res (-) – Residence time in the reactor in s

  • X0 (-) – Composition of the input gas in mole fraction

  • T0_C (-) – Temperature of the input gas in °C

  • P_bar (-) – Pressure in the reactor in bar

  • mechanism (-) – Name of the mechanism file

  • dt_short (-) – Short time step for the simulation

  • conversion_target (-) – Target for the conversion of the input gas. Default is None.

  • H2_yield_target (-) – Target for the H2 yield. Default is None.

  • termination (-) – Termination condition for the simulation. Default is ‘residence_time’. Can be ‘residence_time’ or ‘reactor_length’.

  • L_reactor (-) – Length of the reactor in m. Default is 0. If termination is ‘reactor_length’, L_reactor must be defined.

  • S_reactor (-) – Section of the reactor in m2. Default is 1. If termination is ‘reactor_length’, S_reactor must be defined.

  • name (-) – Name of the reactor. Default is None.

  • verbose (int) – if >0, print infos. If >=2, add calculation progress bar.

Returns:

  • “gas”: ct.Solution

    Gas object with the final state of the reactor

  • ”states”: ct.SolutionArray

    Solution array with the states of the reactor at each time step

Return type:

a dictionary

bloc.reactors.mix_two_streams(gas_1, qm_1, gas_2, qm_2)#

Mix two streams of gases and return the resulting gas object.

bloc.reactors.compute_reactor_length(qm_in, S_in, states)#

qm_in = v * S_in * density <=> v = qm_in / (S_in*density).

qm_in: float

Mass flow rate of the input gas in kg/s

S_in: float

Section of the reactor in m2

states: ct.SolutionArray

contains the results of the kinetic simulation

bloc.reactors.compute_reactor_dimensions_with_form_factor(qm_in, f_factor, states)#
bloc.reactors.get_reactor_mass(L_reactor, d_reactor, e_insul_layers, rho_insul_layers, verbose=False)#

Compute the mass of the reactor based on its dimensions and the insulation layers.

Parameters:
  • L_reactor (float) – Length of the reactor in m.

  • d_reactor (float) – Diameter of the reactor in m.

  • e_insul_layers (list of float) – Thickness of the insulation layers in m.

  • rho_insul_layers (list of float) – Density of the insulation layers in kg/m3.

  • verbose (bool) – If True, print detailed information about the function call.

Returns:

Mass of the reactor in kg.

Return type:

float

bloc.reactors.thermal_resistance_radial_conduction(diameter, length, conductivity, e_insulation)#

Compute the thermal resistance of the radial heat transfer in cylindrical geometry.

bloc.reactors.thermal_resistance_CC(h_CC, S)#

Compute the thermal resistance of the conducto-convective heat transfer.

bloc.reactors.g = 9.81#
bloc.reactors.sigma_SB = 5.67e-08#
bloc.reactors.nu_air = 1.589e-05#
bloc.reactors.alpha_air = 2.25e-05#
bloc.reactors.conductivity_air = 0.0263#
bloc.reactors.Pr_air = 0.707#
bloc.reactors.compute_hCC_natConv_vertCyl(diameter, length, T_wall_C, T_amb_C, conductivity=conductivity_air, nu=nu_air, alpha=alpha_air, Pr=Pr_air)#

Compute the h coefficient for natural convection at the external wall of a vertical cylinder.

bloc.reactors.compute_hrad_linearRadiation(eps, Tmoy)#
bloc.reactors.compute_heat_losses_linear(Tr_C, Tamb_C, R_list)#

Solve the heat transfer problem in the linear approximation.

Tr_C: float

Temperature inside the reactor in °C.

Tamb_C: float

Ambient temperature in °C.

R_list: list

List of the thermal resistances in W/K, from the inside to the outside of the reactor.

S_list: list

List of the surface areas in m2, from the inside to the outside of the reactor.

Returns:

Dictionary containing the heat losses in W and the temperatures at each layer. - “Phi”: float

Heat losses in W

  • ”T_reac”: float

    Temperature inside the reactor in °C

  • ”T_wall_ext”: float

    Temperature of the external wall in °C

  • ”T_amb”: float

    Ambient temperature in °C

Return type:

dict

bloc.reactors.compute_gas2wall_radiative_flux(T_g, T_w, kappa_grey, D)#

Compute the radiative flux from the gas to the wall in W/m2.

Assumptions:
  • temperature and composition are uniform in the reactor

  • reactor is a sphere

  • grey body approximation

  • walls are black (covered by soot) –> no reflection

Reference: Modest, p… TODO.

T_g: float

Gas temperature in K

T_w: float

Wall temperature in K

kappa_grey: float

Absorption coefficient in 1/m in the grey gas approximation

D: float

Reactor diameter in m

bloc.reactors.compute_recovered_power(qm_torch, T_torch_input_C, X_torch_input, qm_2nd_inj, T_2nd_inj_C, X_2nd_inj, P_bar, T_amb, gas_reac, verbose=False)#

Compute the power recovered from the preheating of the torch and the second injection.

bloc.reactors.compute_residual_heat(P_recovered, gas, qm_tot, T_amb_C, verbose=False)#