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

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.

LagrangianPFRReactor

Lagrangian fluid-parcel reactor for tube furnace plug-flow simulation.

LagrangianPFRNetwork

Network that drives a LagrangianPFRReactor through a tube.

Functions#

get_H2_yield(states)

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

get_carbon_yield(gas[, n_C_min])

Compute carbon yield (solid C mass / total C mass in feed).

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).

solve_fixed_position_temperature_profile_pfr(qm_kgs, ...)

Solve kinetics in a PFR with a prescribed temperature profile T(x).

solve_fixed_time_temperature_profile_pfr(...[, dt, ...])

Solve kinetics in a PFR with a prescribed temperature profile T(t).

mix_two_streams(gas_1, qm_1, gas_2, qm_2)

Mix two streams of gases and return the resulting cantera.Quantity 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.

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)

Compute residual heat after recovery in the exchanger.

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 / volumetric_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

class bloc.reactors.LagrangianPFRReactor(gas, *args, clone=False, diameter=0.1, kappa_grey=0.0, mass_flow_rate=0.0, **kwargs)#

Bases: cantera.ExtensibleIdealGasConstPressureMoleReactor

Lagrangian fluid-parcel reactor for tube furnace plug-flow simulation.

This reactor represents a moving fluid parcel advected through the tube. Use LagrangianPFRNetwork to run a simulation — it manages the position-based loop and exposes an advance() interface consistent with cantera.ReactorNet.

The wall contact area is computed dynamically from the current reactor volume as A_internal = (4 / D) * V, which gives the correct surface-to-volume ratio for a circular tube. Because the parcel is a closed system (no mass flow in/out), its volume expands with temperature at constant pressure, exactly as a real fluid element does.

The reactor is intended to be used without MassFlowController connections — just a bare ReactorNet. A Cantera AdaptivePreconditioner can be attached to the network for a significant speedup with large mechanisms (CRECK: ~5–20×).

Parameters:
  • gas (cantera.Solution) – Gas at the initial thermodynamic state.

  • diameter (float) – Tube inner diameter [m].

  • kappa_grey (float, optional) – Grey-gas absorption coefficient [1/m]. Default 0 (no radiation).

  • mass_flow_rate (float, optional) – Mass flow rate [kg/s]. Used for Re / h_conv calculation only.

Notes

T_wall_K and x_position are updated by LagrangianPFRNetwork before each CVODE step. Within one step, after_eval is called several times; all calls reuse the same T_wall_K, which is accurate to first order in position.

See also

LagrangianPFRNetwork

Network that drives this reactor through the tube.

diameter = 0.1#
kappa_grey = 0.0#
mass_flow_rate = 0.0#
T_wall_K: float | None = None#
x_position: float = 0.0#
after_eval(t, LHS, RHS)#

Add wall heat transfer to the temperature equation Right-Hand Side (RHS).

The wall contact area is (4/D) * V where V is the current reactor volume, giving the correct specific surface area for a circular tube. All intermediate calls within one CVODE step use the same T_wall_K (first-order-accurate in position).

After each call the computed values are stored on the reactor for external inspection or diagnostics:

  • _last_h_conv [W/m²/K] — convection coefficient

  • _last_Q_conv [W] — convective heat input to the parcel

  • _last_Q_rad [W] — radiative heat input to the parcel

class bloc.reactors.LagrangianPFRNetwork(reactor, total_length, mass_flow_rate, wall_T_fn, rtol=0.0001, atol=1e-12, use_preconditioner=True)#

Network that drives a LagrangianPFRReactor through a tube.

Wraps a cantera.ReactorNet and runs the position-based advance loop internally. The public interface mirrors cantera.ReactorNet: advance() returns None and results are read from instance attributes afterwards.

Parameters:
  • reactor (LagrangianPFRReactor) – Reactor at the initial (inlet) thermodynamic state, with volume already set to the full tube cross-section times total_length.

  • total_length (float) – Tube length [m].

  • mass_flow_rate (float) – Mass flow rate [kg/s].

  • wall_T_fn (callable) – wall_T_fn(x: float) -> float — wall temperature [K] at position x [m] from the tube inlet.

  • rtol (float, optional) – Relative tolerance for the CVODE integrator. Default 1e-4.

  • atol (float, optional) – Absolute tolerance for the CVODE integrator. Default 1e-12.

  • use_preconditioner (bool, optional) – Attach an cantera.AdaptivePreconditioner to the network. Default True.

net#

Underlying reactor network. Configure it (e.g. set a custom preconditioner, change tolerances, or override max_time_step) before calling advance(). A conservative max_time_step is set automatically so that CVODE cannot skip the thermal entry length in a single step (important for non-reactive carrier gases such as N2 where there is no chemical stiffness to naturally limit the step size).

Type:

cantera.ReactorNet

n_steps#

Number of CVODE steps taken during the last advance() call.

Type:

int

E_conv#

Total convective energy transferred from wall to gas [J].

Type:

float

E_rad#

Total radiative energy transferred from wall to gas [J].

Type:

float

x_arr#

Parcel positions [m] at each step. Set only when record_profiles is True in advance(), otherwise None.

Type:

numpy.ndarray or None

T_wall_arr#

Wall temperature [K] at each recorded position.

Type:

numpy.ndarray or None

T_gas_arr#

Gas temperature [K] at each recorded position.

Type:

numpy.ndarray or None

Y_Cs_arr#

Solid carbon mass fraction at each recorded position.

Type:

numpy.ndarray or None

Notes

This class uses a custom Lagrangian advance loop rather than a standard ReactorNet simulation. An alternative architecture — building a Wall-based steady-state network solved directly with ReactorNet — is possible but not implemented here; see build_visualization_network() for details.

Examples

>>> reactor = LagrangianPFRReactor(
...     gas, diameter=0.1, kappa_grey=0.5, mass_flow_rate=mdot
... )
>>> reactor.volume = np.pi / 4 * 0.1**2 * total_length
>>> network = LagrangianPFRNetwork(reactor, total_length, mdot, wall_T_fn)
>>> network.advance()
>>> T_outlet = reactor.phase.T
reactor#
total_length#
mass_flow_rate#
wall_T_fn#
net#
n_steps: int = 0#
E_conv: float = 0.0#
E_rad: float = 0.0#
h_conv_avg: float#
x_arr: numpy.ndarray | None = None#
t_arr: numpy.ndarray | None = None#
T_wall_arr: numpy.ndarray | None = None#
T_gas_arr: numpy.ndarray | None = None#
Y_Cs_arr: numpy.ndarray | None = None#
X_arr: numpy.ndarray | None = None#
species_names: list[str] | None = None#
Fo_D_arr: numpy.ndarray | None = None#
Fo_D_min: float#
h_conv_arr: numpy.ndarray | None = None#
alpha_arr: numpy.ndarray | None = None#
u_arr: numpy.ndarray | None = None#
k_arr: numpy.ndarray | None = None#
Re_arr: numpy.ndarray | None = None#
Pr_arr: numpy.ndarray | None = None#
advance(verbose=False, dx_print=0.05, record_profiles=False)#

Advance the parcel from the tube inlet to the outlet.

At each CVODE step the parcel position is updated from the local gas density and the mass flow rate, and reactor.T_wall_K / reactor.x_position are set accordingly. The loop ends when the parcel has travelled total_length.

Returns None — consistent with cantera.ReactorNet.advance(). Read results from instance attributes after the call: n_steps, E_conv, E_rad, and (when record_profiles=True) x_arr, T_wall_arr, T_gas_arr.

Parameters:
  • verbose (bool, optional) – Print heat-transfer diagnostics every dx_print m. Default False.

  • dx_print (float, optional) – Spatial interval between diagnostic lines [m]. Default 0.05 m.

  • record_profiles (bool, optional) – When True, record (x, T_wall, T_gas, Y_Cs) at every CVODE step and store as numpy arrays in x_arr, T_wall_arr, T_gas_arr, Y_Cs_arr. Default False.

build_visualization_network(P_total_W=None)#

Build a visualization-only ReactorNet representing the tube topology.

Creates Inlet [MFC] Tube [MFC] Outlet with T_wall [Wall] Tube. This network is never solved; it exists solely for diagram generation via ct.ReactorNet.draw().

When P_total_W is provided (convective + radiative power from the Lagrangian simulation), the Wall heat flux is set so that Cantera’s draw() reports the actual delivered power instead of = 0 W.

Notes

An alternative simulation architecture — building a proper Wall-based steady-state network solved directly with ReactorNet — is possible but not implemented here. In that approach a single IdealGasConstPressureReactor would be connected to a Reservoir at T_wall via a Wall with a prescribed heat transfer coefficient derived from the Lagrangian h_conv. Comparing the two approaches provides a useful cross-validation of the Lagrangian result.

Parameters:

P_total_W (float or None) – Total wall-to-gas power [W] = (E_conv + E_rad) / t_res. When given, the Wall heat flux is set to P_total_W / A_lateral so the diagram shows the actual delivered power.

Returns:

Initialised (but unsolved) network ready for draw().

Return type:

cantera.ReactorNet

NETWORK_DIAGRAM_NOTE: str = 'Node temperatures in this diagram reflect the outlet gas state of the Lagrangian simulation...#
draw(title='Tube Furnace Wall Temperature Profile', n_points=200)#

Return a schematic matplotlib Figure of the tube furnace geometry.

The upper panel shows a colour-coded tube rectangle (cold→hot→cold). The lower panel shows the wall temperature profile (and the gas temperature profile when advance() was already called with record_profiles=True).

This method does not require advance() to have been called; it relies only on wall_T_fn and total_length.

Return type:

matplotlib.figure.Figure

plot_temperature_profile(title='Temperature Profile Tube Furnace')#

Plot gas and wall temperature vs position with residence-time axis.

Requires advance() to have been called with record_profiles=True.

Return type:

matplotlib.figure.Figure

plot_species_profiles(species=None, title='Species Profile Tube Furnace')#

Plot mole-fraction profiles vs position with residence-time axis.

Requires advance() to have been called with record_profiles=True.

Parameters:

species (list of str or None) – Species to plot. When None, the six species with the highest peak mole fraction (≥ 0.1 %) are selected automatically.

Return type:

matplotlib.figure.Figure

bloc.reactors.get_H2_yield(states)#

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

Parameters:
  • states (ct.SolutionArray) – Solution array with the states of the reactor at each time step. states[0] and states[-1] must be the initial and final states of the reactor. Mass is assumed to be conserved between states[0] and states[-1].

  • as (The H2 yield is defined)

  • subtracted (The initial amount of H2 in the feedstock must be)

  • produced. (because it is not)

  • Thus (H2_yield = (Y_H2_final - Y_H2_initial) / (Y_H_tot - Y_H2_initial))

  • math:: (..) – H2_yield = \frac{Y_{H2,final} - Y_{H2,initial}}{Y_{H,total} - Y_{H2,initial}}

bloc.reactors.get_carbon_yield(gas, n_C_min=300)#

Compute carbon yield (solid C mass / total C mass in feed).

Assumes no solid carbon in the input gas. Carbon yield = mass fraction of solid carbon / total mass fraction of carbon element in the input gas.

Parameters:
  • gas (cantera.Solution or cantera.SolutionArray) – Gas object with composition set.

  • n_C_min (int, optional) – Minimum number of carbon atoms to treat a species as solid carbon. Default 300. Use 24 for mechanisms that represent soot as large PAH only.

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, P_in_kW_avg=0.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 (-) – Instantaneous power input in kW. Applied once at t = 0 as an inlet specific-enthalpy pulse P_in_kW * 1e3 / qm_kgs (J/kg), before time integration starts.

  • 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.

  • P_in_kW_avg (-) – Average power input/removal in kW distributed uniformly over t_res. At each integration step dt, the specific enthalpy is adjusted by dh = P_in_kW_avg * 1e3 / qm_kgs * dt / t_res (J/kg). Positive values add heat; negative values remove heat.

  • 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.solve_fixed_position_temperature_profile_pfr(qm_kgs, inlet_composition, P_bar, mechanism, temperature_profile_position, dx=0.01, L_reactor=0.0, d_reactor=1.0, name=None, verbose=2)#

Solve kinetics in a PFR with a prescribed temperature profile T(x).

Integrates along the reactor axis with fixed spatial step dx. At each step, the reactor advances by dt = dx / velocity (from mass flow and cross-section), then temperature is set from the profile. No energy equation; T is imposed.

Parameters:
  • qm_kgs (float) – Mass flow rate (kg/s).

  • inlet_composition (str) – Inlet composition (Cantera format, mole fractions).

  • P_bar (float) – Pressure (bar).

  • mechanism (str, optional) – Mechanism name or path; must not be None.

  • temperature_profile_position (np.ndarray) – Shape (2, n): row 0 = axial positions (m), row 1 = temperature (K).

  • dx (float, optional) – Spatial integration step (m). Default 0.01.

  • L_reactor (float, optional) – Reactor length (m); used for reporting. Default 0.0.

  • d_reactor (float, optional) – Reactor diameter (m); used for velocity and cross-section. Default 1.0.

  • name (str, optional) – Reactor name for logging.

  • verbose (int, optional) – Logging level. Default 2.

Returns:

“gas”: Cantera Solution; “states”: SolutionArray with extra “t”, “x”. If mechanism is CRECK_Nobili2024.yaml, also “reaction_rates_by_class” and “mass_carbon_rates_by_class”.

Return type:

dict

bloc.reactors.solve_fixed_time_temperature_profile_pfr(inlet_composition, P_bar, mechanism, temperature_profile_time, dt=1e-06, name=None, verbose=2)#

Solve kinetics in a PFR with a prescribed temperature profile T(t).

Integrates in time with fixed step dt. At each step the reactor advances by dt, then temperature is set from the profile. No energy equation; T is imposed. No reactor geometry or velocity; purely time-residence integration.

Parameters:
  • inlet_composition (str) – Inlet composition (Cantera format, mole fractions).

  • P_bar (float) – Pressure (bar).

  • mechanism (str, optional) – Mechanism name or path; must not be None.

  • temperature_profile_time (np.ndarray) – Shape (2, n): row 0 = residence times (s), row 1 = temperature (K).

  • dt (float, optional) – Time integration step (s). Default 1e-6.

  • name (str, optional) – Reactor name for logging.

  • verbose (int, optional) – Logging level. Default 2.

Returns:

“gas”: Cantera Solution; “states”: SolutionArray with extra “t”. If mechanism is CRECK_Nobili2024.yaml, also “reaction_rates_by_class” and “mass_carbon_rates_by_class”.

Return type:

dict

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

Mix two streams of gases and return the resulting cantera.Quantity 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.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.

\[P_{\mathrm{recovered}} = P_{\mathrm{preheat,torch}} + P_{\mathrm{preheat,2nd}} P_{\mathrm{preheat}} = \dot{m} \, (h_{\mathrm{in}} - h_0)\]
bloc.reactors.compute_residual_heat(P_recovered, gas, qm_tot, T_amb_C, verbose=False)#

Compute residual heat after recovery in the exchanger.

\[P_{\mathrm{heat,tot}} = \dot{m}_{\mathrm{tot}} \, (h_{\mathrm{reactor\,out}} - h_{\mathrm{cold\,out}}) P_{\mathrm{residual}} = P_{\mathrm{heat,tot}} - P_{\mathrm{recovered}}\]