Fixed Bed 1D Reactor#

The IDAES Fixed Bed 1D Reactor (FixedBed1D) model represents a unit operation where a gas stream passes through a solid phase bed in a linear reactor vessel. The FixedBed1D mathematical model is a 1-D time variant model with two phases (gas and solid). The model captures the gas-solid interaction between both phases through reaction/adsorption, mass and heat transfer.

Assumptions:

  • The radial concentration and temperature gradients are assumed to be negligible.

  • The reactor is assumed to be adiabatic.

Requirements:

  • Property package contains temperature and pressure variables.

Degrees of Freedom#

FixedBed1D Reactors generally have at least 2 (or more) degrees of freedom, consisting of design and operating variables. The design variables of reactor length and diameter are typically the minimum variables to be fixed.

Model Structure#

The core FixedBed1D unit model consists of one ControlVolume1DBlock Block named gas_phase, which has one Inlet Port (named gas_inlet) and one Outlet Port (named gas_outlet). As there are no flow variables in the solid_phase a ControlVolume1DBlock Block is not used, instead indexed state blocks are used to access the solid thermophysical properties whilst the mass and energy balance equations of the solid are written in the core FixedBed1D unit model.

Construction Arguments#

The IDAES FixedBed1D model has construction arguments specific to the whole unit and to the individual regions.

Arguments that are applicable to the FixedBed1D unit as a whole are:

  • finite_elements - sets the number of finite elements to use when discretizing the spatial domains (default = 10).

  • length_domain_set - sets the list of point to use to initialize a new ContinuousSet (default = [0.0, 1.0]).

  • transformation_method - sets the discretization method to use by the Pyomo TransformationFactory to transform the spatial domain (default = dae.finite_difference):

    • dae.finite_difference - finite difference method.

    • dae.collocation - orthogonal collocation method.

  • transformation_scheme - sets the scheme to use when transforming a domain. Selected schemes should be compatible with the transformation_method chosen (default = None):

    • None - defaults to “BACKWARD” for finite difference transformation method and to “LAGRANGE-RADAU” for collocation transformation method

    • BACKWARD - use a finite difference transformation method.

    • FORWARD - use a finite difference transformation method.

    • LAGRANGE-RADAU - use a collocation transformation method.

  • collocation_points - sets the number of collocation points to use when discretizing the spatial domains (default = 3, collocation methods only).

  • flow_type - indicates the flow arrangement within the unit to be modeled. Options are:

    • ‘forward_flow’ - (default) gas flows in the forward direction (from x=0 to x=1)

    • ‘reverse_flow’ - gas flows in the reverse direction (from x=1 to x=0).

  • material_balance_type - indicates what type of energy balance should be constructed (default = MaterialBalanceType.componentTotal).

    • MaterialBalanceType.componentTotal - use total component balances.

    • MaterialBalanceType.total - use total material balance.

  • energy_balance_type - indicates what type of energy balance should be constructed (default = EnergyBalanceType.enthalpyTotal).

    • EnergyBalanceType.none - excludes energy balances.

    • EnergyBalanceType.enthalpyTotal - single enthalpy balance for material.

  • momentum_balance_type - indicates what type of momentum balance should be constructed (default = MomentumBalanceType.pressureTotal).

    • MomentumBalanceType.none - exclude momentum balances.

    • MomentumBalanceType.pressureTotal - single pressure balance for material.

  • has_pressure_change - indicates whether terms for pressure change should be constructed (default = True).

    • True - include pressure change terms.

    • False - exclude pressure change terms.

  • pressure_drop_type - indicates what type of pressure drop correlation should be used (default = “ergun_correlation”).

    • “ergun_correlation” - use the Ergun equation.

    • “simple_correlation” - use a simplified pressure drop correlation.

Arguments that are applicable to the gas phase:

  • property_package - property package to use when constructing gas phase Property Blocks (default = ‘use_parent_value’). This is provided as a Physical Parameter Block by the Flowsheet when creating the model. If a value is not provided, the ControlVolume Block will try to use the default property package if one is defined.

  • property_package_args - set of arguments to be passed to the gas phase Property Blocks when they are created (default = ‘use_parent_value’).

  • reaction_package - reaction package to use when constructing gas phase Reaction Blocks (default = None). This is provided as a Reaction Parameter Block by the Flowsheet when creating the model. If a value is not provided, the ControlVolume Block will try to use the default property package if one is defined.

  • reaction_package_args - set of arguments to be passed to the gas phase Reaction Blocks when they are created (default = None).

  • has_equilibrium_reactions - sets flag to indicate if terms of equilibrium controlled reactions should be constructed (default = False).

Arguments that are applicable to the solid phase:

  • property_package - property package to use when constructing solid phase Property Blocks (default = ‘use_parent_value’). This is provided as a Physical Parameter Block by the Flowsheet when creating the model. If a value is not provided, the Indexed State Blocks will try to use the default property package if one is defined.

  • property_package_args - set of arguments to be passed to the solid phase Property Blocks when they are created (default = ‘use_parent_value’).

  • reaction_package - reaction package to use when constructing solid phase Reaction Blocks (default = None). This is provided as a Reaction Parameter Block by the Flowsheet when creating the model. If a value is not provided, the Indexed State Blocks will try to use the default property package if one is defined.

  • reaction_package_args - set of arguments to be passed to the solid phase Reaction Blocks when they are created (default = None).

  • has_equilibrium_reactions - sets flag to indicate if terms of equilibrium controlled reactions should be constructed (default = False).

Additionally, FixedBed1D units have the following construction arguments which are passed to all the ControlVolume1DBlock and state Blocks and are always specified to True.

Argument

Value

dynamic

True

has_holdup

True

Constraints#

In the following, the subscripts \(g\) and \(s\) refer to the gas and solid phases, respectively. In addition to the constraints written by the control_volume Block, FixedBed1D units write the following Constraints:

Geometry Constraints#

Area of the reactor bed:

\[A_{bed} = \pi \left( \frac{ D_{bed} }{ 2 } \right)^2\]

Area of the gas domain:

\[A_{g,t,x} = \varepsilon A_{bed}\]

Area of the solid domain:

\[A_{s,t,x} = (1 - \varepsilon) A_{bed}\]

Hydrodynamic Constraints#

Superficial velocity of the gas:

\[u_{g,t,x} = \frac{ F_{mol,g,t,x} }{ A_{bed} \rho_{mol,g,t,x} }\]

Pressure drop:

The constraints written by the FixedBed1D model to compute the pressure drop (if has_pressure_change is ‘True’) in the reactor depend upon the construction arguments chosen:

If pressure_drop_type is simple_correlation:

\[- \frac{ dP_{g,t,x} }{ dx } = 0.2 \left( \rho_{mass,s,t,x} - \rho_{mass,g,t,x} \right) u_{g,t,x}\]

If pressure_drop_type is ergun_correlation:

\[- \frac{ dP_{g,t,x} }{ dx } = \frac{ 150 \mu_{g,t,x} {\left( 1 - \varepsilon \right)}^{2} u_{g,t,x} }{ \varepsilon^{3} d_{p}^2 } + \frac{ 1.75 \left( 1 - \varepsilon \right) \rho_{mass,g,t,x} u_{g,t,x}^{2} }{ \varepsilon^{3} d_{p} }\]

Reaction Constraints#

If gas_phase_config.reaction_package is not ‘None’:

Gas phase reaction extent:

\[\xi_{g,t,x,r} = r_{g,t,x,r} A_{g,t,x}\]

If solid_phase_config.reaction_package is not ‘None’:

Gas phase heterogeneous rate generation/consumption:

\[M_{g,t,x,p,j} = A_{s,t,x} \sum_{r}^{reactions} {\nu_{s,j,r} r_{s,t,x,r}}\]

Dimensionless numbers, mass and heat transfer coefficients#

Sum of component mass fractions in solid phase:

\[{\Sigma_{j}} {x_{mass,s,t,x,j}} = 1 \space \space \space \forall \space t \forall \space x\]

Solid material holdup:

\[J_{mass,s,t,x,j} = A_{s,t,x} \rho_{mass,s,t,x} x_{mass,s,t,x,j}\]

If solid_phase_config.reaction_package is not ‘None’:

Solid material accumulation:

\[{\dot{J}}_{mass,s,t,x,j} = A_{s,t,x} {MW}_{s,j} {\Sigma}_{r} {r_{s,t,x,r} \nu_{s,j,r}}\]

If solid_phase_config.reaction_package is ‘None’:

Solid material accumulation:

\[{\dot{J}}_{mass,s,t,x,j} = 0\]

If energy_balance_type not EnergyBalanceType.none:

Solid energy holdup:

\[q_{energy,s,t,x} = A_{s,t,x} \rho_{mass,s,t,x} H_{mass,s,t,x}\]

Solid energy accumulation:

\[{\dot{q}}_{energy,s,t,x} = H_{s,t,x} - A_{s,t,x} {\Sigma}_{r} {r_{s,t,x,r} H_{rxn,s,t,x,r}}\]

If energy_balance_type is EnergyBalanceType.none:

Isothermal solid phase:

\[T_{s,t,x} = T_{s,t=0,x}\]

Isothermal gas phase:

\[T_{g,t,x} = T_{s,t=0,x}\]

If energy_balance_type not EnergyBalanceType.none:

Particle Reynolds number:

\[Re_{p,t,x} = \frac{ u_{g,t,x} \rho_{mass,g,t,x} }{ \mu_{g,t,x} d_{p}}\]

Prandtl number:

\[Pr_{t,x} = \frac{ c_{p,t,x} \mu_{g,t,x} }{ k_{g,t,x} }\]

Particle Nusselt number:

\[Nu_{p,t,x} = 2 + 1.1 Pr_{t,x}^{1/3} \left| Re_{p,t,x} \right|^{0.6}\]

Particle to fluid heat transfer coefficient

\[h_{gs,t,x} d_{p} = Nu_{p,t,x} k_{g,t,x}\]

Gas phase - gas to solid heat transfer:

\[H_{g,t,x} = - \frac{ 6 } { d_{p} } h_{gs,t,x} \left( T_{g,t,x} - T_{s,t,x} \right) A_{s,t,x}\]

Solid phase - gas to solid heat transfer:

\[H_{s,t,x} = \frac{ 6 } { d_{p} } h_{gs,t,x} \left( T_{g,t,x} - T_{s,t,x} \right) A_{s,t,x}\]

List of Variables#

Variable

Description

Reference to

\(A_{bed}\)

Reactor bed cross-sectional area

bed_area

\(A_{g,t,x}\)

Gas phase area (interstitial cross-sectional area)

gas_phase.area

\(A_{s,t,x}\)

Solid phase area

solid_phase.area

\(c_{p,t,x}\)

Gas phase heat capacity (constant \(P\))

gas_phase.properties.cp_mass

\(D_{bed}\)

Reactor bed diameter

bed_diameter

\(F_{mol,g,t,x}\)

Total molar flow rate of gas

gas_phase.properties.flow_mol

\(H_{g,t,x}\)

Gas to solid heat transfer term, gas phase

gas_phase.heat

\(H_{s,t,x}\)

Gas to solid heat transfer term, solid phase

solid_phase.heat

\(H_{rxn,s,t,x,r}\)

Solid phase heat of reaction

solids.reactions.dh_rxn

\(h_{gs,t,x}\)

Gas-solid heat transfer coefficient

gas_solid_htc

\(J_{mass,s,t,x,j}\)

Material holdup, solid phase

solids.solids_material_holdup

\({\dot{J}}_{mass,s,t,x,j}\)

Material accumulation, solid phase

solids.solids_material_accumulation

\(k_{g,t,x}\)

Gas thermal conductivity

gas_phase.properties.therm_cond

\(L_{bed}\)

Reactor bed height

bed_height

\(M_{g,t,x,p,j}\)

Rate generation/consumption term, gas phase

gas_phase.mass_transfer_term

\(Nu_{p,t,x}\)

Particle Nusselt number

Nu_particle

\(dP_{g,t,x}\)

Total pressure derivative w.r.t. \(x\) (axial position)

gas_phase.deltaP

\(Pr_{t,x}\)

Prandtl number

Pr

\(q_{energy,s,t,x}\)

Energy holdup, solid phase

solids.energy_material_holdup

\({\dot{q}}_{energy,s,t,x}\)

Energy accumulation, solid phase

solids.solids_energy_accumulation

\(r_{g,t,x,r}\)

Gas phase reaction rate

gas_phase.reactions.reaction_rate

\(r_{s,t,x,r}\)

Solid phase reaction rate

solid_phase.reactions.reaction_rate

\(Re_{p,t,x}\)

Particle Reynolds number

Re_particle

\(T_{g,t,x}\)

Gas phase temperature

gas_phase.properties.temperature

\(T_{s,t,x}\)

Solid phase temperature

solid_phase.properties.temperature

\(u_{g,t,x}\)

Superficial velocity of the gas

velocity_superficial_gas

Greek letters

\(\varepsilon\)

Reactor bed voidage

bed_voidage

\(\mu_{g,t,x}\)

Dynamic viscosity of gas mixture

gas_phase.properties.visc_d

\(\xi_{g,t,x,r}\)

Gas phase reaction extent

gas_phase.rate_reaction_extent

\(\rho_{mass,g,t,inlet}\)

Density of gas mixture

gas_phase.properties.dens_mass

\(\rho_{mass,s,t,inlet}\)

Density of solid particles

solid_phase.properties.dens_mass_particle

\(\rho_{mol,g,t,x}\)

Molar density of the gas

gas_phase.properties.dens_mole

List of Parameters#

Parameter

Description

Reference to

\(d_{p}\)

Solid particle diameter

solid_phase.properties._params.particle_dia

\(\nu_{s,j,r}\)

Stoichiometric coefficients

solid_phase.reactions.rate_reaction_stoichiometry

Initialization#

The initialization method for this model will save the current state of the model before commencing initialization and reloads it afterwards. The state of the model will be the same after initialization, only the initial guesses for unfixed variables will be changed.

The model allows for the passing of a dictionary of values of the state variables of the gas and solid phases that can be used as initial guesses for the state variables throughout the time and spatial domains of the model. This is optional but recommended. A typical guess could be values of the gas and solid inlet port variables at time \(t=0\).

The model initialization proceeds through a sequential hierarchical method where the model equations are deactivated at the start of the initialization routine, and the complexity of the model is built up through activation and solution of various sub-model blocks and equations at each initialization step. At each step the model variables are updated to better guesses obtained from the model solution at that step.

The initialization routine proceeds as follows:

  • Step 1: Initialize the thermo-physical and transport properties model blocks.

  • Step 2: Initialize the hydrodynamic properties.

  • Step 3a: Initialize mass balances without reactions and pressure drop.

  • Step 3b: Initialize mass balances with reactions and without pressure drop.

  • Step 3c: Initialize mass balances with reactions and pressure drop.

  • Step 4: Initialize energy balances.

Block Triangularization Initialization#

The Block Triangularization initialization routine initializes the 1DFixedBed model in the following steps:

  • Step 1: Convert the system of equations from a partial differential set of equations (PDAE) to an algebraic set of equations (AE) by deactivating

    the discretization equations and sum_component_equations (if present), and fixing the state variables to a guess (typically inlet conditions).

  • Step 2: Decompose the AE into strongly connected components via block triangularization and solve individual blocks.

  • Step 3: Revert the fixed variables and deactivated constraints to their original states.

The routine allows for the passing of a dictionary of values of the state variables of the gas and solid phases that can be used as initial guesses for the state variables throughout the time and spatial domains of the model. This is optional but recommended. A typical guess could be values of the gas and solid state variables at time \(t=0\).

FixedBed1D Class#

class idaes.models_extra.gas_solid_contactors.unit_models.fixed_bed_1D.FixedBed1D(*args, **kwds)#
Parameters:
  • rule (function) – A rule function or None. Default rule calls build().

  • concrete (bool) – If True, make this a toplevel model. Default - False.

  • ctype (class) –

    Pyomo ctype of the block. Default - pyomo.environ.Block

    Config args

    dynamic

    Indicates whether this model will be dynamic or not, default = useDefault. Valid values: { useDefault - get flag from parent (default = False), True - set as a dynamic model, False - set as a steady-state model.}

    has_holdup

    Indicates whether holdup terms should be constructed or not. Must be True if dynamic = True, default - False. Valid values: { useDefault - get flag from parent (default = False), True - construct holdup terms, False - do not construct holdup terms}

    finite_elements

    Number of finite elements to use when discretizing length domain (default=20)

    length_domain_set

    length_domain_set - (optional) list of point to use to initialize a new ContinuousSet if length_domain is not provided (default = [0.0, 1.0])

    transformation_method

    Method to use to transform domain. Must be a method recognized by the Pyomo TransformationFactory, default - “dae.finite_difference”. Valid values: { “dae.finite_difference” - Use a finite difference transformation method, “dae.collocation” - use a collocation transformation method}

    transformation_scheme

    Scheme to use when transforming domain. See Pyomo documentation for supported schemes, default - None. Valid values: { None - defaults to “BACKWARD” for finite difference transformation method, and to “LAGRANGE- RADAU” for collocation transformation method, “BACKWARD” - Use a finite difference transformation method, “FORWARD”” - use a finite difference transformation method, “LAGRANGE-RADAU”” - use a collocation transformation method}

    collocation_points

    Number of collocation points to use per finite element when discretizing length domain (default=3)

    flow_type

    Flow configuration of Fixed Bed default - “forward_flow”. Valid values: { “forward_flow” - gas flows from 0 to 1, “reverse_flow” - gas flows from 1 to 0.}

    material_balance_type

    Indicates what type of mass balance should be constructed, default - MaterialBalanceType.componentTotal. Valid values: { MaterialBalanceType.none - exclude material balances, MaterialBalanceType.componentPhase - use phase component balances, MaterialBalanceType.componentTotal - use total component balances, MaterialBalanceType.elementTotal - use total element balances, MaterialBalanceType.total - use total material balance.}

    energy_balance_type

    Indicates what type of energy balance should be constructed, default - EnergyBalanceType.enthalpyTotal. Valid values: { EnergyBalanceType.none - exclude energy balances, EnergyBalanceType.enthalpyTotal - single enthalpy balance for material, EnergyBalanceType.enthalpyPhase - enthalpy balances for each phase, EnergyBalanceType.energyTotal - single energy balance for material, EnergyBalanceType.energyPhase - energy balances for each phase.}

    momentum_balance_type

    Indicates what type of momentum balance should be constructed, default - MomentumBalanceType.pressureTotal. Valid values: { MomentumBalanceType.none - exclude momentum balances, MomentumBalanceType.pressureTotal - single pressure balance for material, MomentumBalanceType.pressurePhase - pressure balances for each phase, MomentumBalanceType.momentumTotal - single momentum balance for material, MomentumBalanceType.momentumPhase - momentum balances for each phase.}

    has_pressure_change

    Indicates whether terms for pressure change should be constructed, default - False. Valid values: { True - include pressure change terms, False - exclude pressure change terms.}

    pressure_drop_type

    Indicates what type of pressure drop correlation should be used, default - “ergun_correlation”. Valid values: { “ergun_correlation” - Use the Ergun equation, “simple_correlation” - Use a simplified pressure drop correlation.}

    gas_phase_config

    gas phase config arguments

    gas_phase_config
    dynamic

    Indicates whether this model will be dynamic or not, default = useDefault. Valid values: { useDefault - get flag from parent (default = False), True - set as a dynamic model, False - set as a steady-state model.}

    has_holdup

    Indicates whether holdup terms should be constructed or not. Must be True if dynamic = True, default - False. Valid values: { useDefault - get flag from parent (default = False), True - construct holdup terms, False - do not construct holdup terms}

    property_package

    Property parameter object used to define property calculations (default = ‘use_parent_value’) - ‘use_parent_value’ - get package from parent (default = None) - a ParameterBlock object

    property_package_args

    A dict of arguments to be passed to the PropertyBlockData and used when constructing these (default = ‘use_parent_value’) - ‘use_parent_value’ - get package from parent (default = None) - a dict (see property package for documentation)

    reaction_package

    Reaction parameter object used to define reaction calculations, default - None. Valid values: { None - no reaction package, ReactionParameterBlock - a ReactionParameterBlock object.}

    reaction_package_args

    A ConfigBlock with arguments to be passed to a reaction block(s) and used when constructing these, default - None. Valid values: { see reaction package for documentation.}

    has_equilibrium_reactions

    Indicates whether terms for equilibrium controlled reactions should be constructed, default - True. Valid values: { True - include equilibrium reaction terms, False - exclude equilibrium reaction terms.}

    solid_phase_config

    solid phase config arguments

    solid_phase_config
    dynamic

    Indicates whether this model will be dynamic or not, default = useDefault. Valid values: { useDefault - get flag from parent (default = False), True - set as a dynamic model, False - set as a steady-state model.}

    has_holdup

    Indicates whether holdup terms should be constructed or not. Must be True if dynamic = True, default - False. Valid values: { useDefault - get flag from parent (default = False), True - construct holdup terms, False - do not construct holdup terms}

    property_package

    Property parameter object used to define property calculations (default = ‘use_parent_value’) - ‘use_parent_value’ - get package from parent (default = None) - a ParameterBlock object

    property_package_args

    A dict of arguments to be passed to the PropertyBlockData and used when constructing these (default = ‘use_parent_value’) - ‘use_parent_value’ - get package from parent (default = None) - a dict (see property package for documentation)

    reaction_package

    Reaction parameter object used to define reaction calculations, default - None. Valid values: { None - no reaction package, ReactionParameterBlock - a ReactionParameterBlock object.}

    reaction_package_args

    A ConfigBlock with arguments to be passed to a reaction block(s) and used when constructing these, default - None. Valid values: { see reaction package for documentation.}

    has_equilibrium_reactions

    Indicates whether terms for equilibrium controlled reactions should be constructed, default - True. Valid values: { True - include equilibrium reaction terms, False - exclude equilibrium reaction terms.}

  • initialize (dict) – ProcessBlockData config for individual elements. Keys are BlockData indexes and values are dictionaries with config arguments as keys.

  • idx_map (function) – Function to take the index of a BlockData element and return the index in the initialize dict from which to read arguments. This can be provided to override the default behavior of matching the BlockData index exactly to the index in initialize.

Returns:

(FixedBed1D) New instance

FixedBed1DData Class#

class idaes.models_extra.gas_solid_contactors.unit_models.fixed_bed_1D.FixedBed1DData(component)[source]#

1D Fixed Bed Reactor Model Class

block_triangularization_initialize(gas_phase_state_args=None, solid_phase_state_args=None, outlvl=0, solver=None, calc_var_kwds=None)[source]#

Block triangularization (BT) initialization routine for 1D FixedBed unit.

The BT initialization routine initializes the unit model in the following steps:

Step 1: Convert the system of equations from a partial differential set of

equations (PDAE) to an algebraic set of equations (AE) by deactivating the discretization equations and sum_component_equations (if present), and fixing the state variables to a guess (typically inlet conditions).

Step 2: Decompose the AE into strongly connected components

via block triangularization and solve individual blocks

Step 3: revert the fixed variables and deactivated constraints to their

original states.

The initialization procedure is wrapped in the TemporarySubsystemManager utility to ensure that the deactivation of constraints and fixing of variables are only temporary and the final model structure is unchanged after initialization.

Keyword Arguments:
  • gas_phase_state_args – a dict of arguments to be passed to the property package(s) to provide an initial state for initialization (see documentation of the specific property package) (default = None).

  • solid_phase_state_args – a dict of arguments to be passed to the property package(s) to provide an initial state for initialization (see documentation of the specific property package) (default = None).

  • outlvl – sets output level of initialization routine

  • optarg – solver options dictionary object (default=None, use default solver options)

  • solver – str indicating which solver to use during initialization (default = None, use default solver)

  • calc_var_kwds – Dictionary of keyword arguments for calculate_variable_from_constraint

Returns:

None

build()[source]#

Begin building model (pre-DAE transformation).

Parameters:

None

Returns:

None

initialize_build(gas_phase_state_args=None, solid_phase_state_args=None, outlvl=0, solver=None, optarg=None)[source]#

Initialization routine for 1DFixedBed unit.

Keyword Arguments:
  • gas_phase_state_args – a dict of arguments to be passed to the property package(s) to provide an initial state for initialization (see documentation of the specific property package) (default = None).

  • solid_phase_state_args – a dict of arguments to be passed to the property package(s) to provide an initial state for initialization (see documentation of the specific property package) (default = None).

  • outlvl – sets output level of initialization routine

  • optarg – solver options dictionary object (default=None, use default solver options)

  • solver – str indicating which solver to use during initialization (default = None, use default solver)

Returns:

None