Hide code cell content
###############################################################################
# The Institute for the Design of Advanced Energy Systems Integrated Platform
# Framework (idaes IP) was produced under the DOE Institute for the
# Design of Advanced Energy Systems (IDAES).
#
# Copyright (c) 2018-2023 by the software owners: The Regents of the
# University of California, through Lawrence Berkeley National Laboratory,
# National Technology & Engineering Solutions of Sandia, LLC, Carnegie Mellon
# University, West Virginia University Research Corporation, et al.
# All rights reserved.  Please see the files COPYRIGHT.md and LICENSE.md
# for full copyright and license information.
###############################################################################

HDA Flowsheet Simulation and Optimization#

Author: Jaffer Ghouse
Maintainer: Tanner Polley
Updated: 2026-1-12

Learning outcomes#

  • Construct a steady-state flowsheet using the IDAES unit model library

  • Connecting unit models in a flowsheet using Arcs

  • Using the SequentialDecomposition tool to initialize a flowsheet with recycle

  • Formulate and solve an optimization problem

    • Defining an objective function

    • Setting variable bounds

    • Adding additional constraints

The general workflow of setting up an IDAES flowsheet is the following:

     1 Importing Modules
     2 Building a Model
     3 Scaling the Model
     4 Specifying the Model
     5 Initializing the Model
     6 Solving the Model
     7 Analyzing and Visualizing the Results
     8 Optimizing the Model

We will complete each of these steps as well as demonstrate analyses on this model through some examples and exercises.

Problem Statement#

Hydrodealkylation is a chemical reaction that often involves reacting an aromatic hydrocarbon in the presence of hydrogen gas to form a simpler aromatic hydrocarbon devoid of functional groups. In this example, toluene will be reacted with hydrogen gas at high temperatures to form benzene via the following reaction:

C6H5CH3 + H2 → C6H6 + CH4

This reaction is often accompanied by an equilibrium side reaction which forms diphenyl, which we will neglect for this example.

This example is based on the 1967 AIChE Student Contest problem as present by Douglas, J.M., Chemical Design of Chemical Processes, 1988, McGraw-Hill.

The flowsheet that we will be using for this module is shown below with the stream conditions. We will be processing toluene and hydrogen to produce at least 370 TPY of benzene. As shown in the flowsheet, there are two flash tanks, F101 to separate out the non-condensibles and F102 to further separate the benzene-toluene mixture to improve the benzene purity. Note that typically a distillation column is required to obtain high purity benzene but that is beyond the scope of this workshop. The non-condensibles separated out in F101 will be partially recycled back to M101 and the rest will be either purged or combusted for power generation.We will assume ideal gas for this flowsheet. The properties required for this module are available in the same directory:

  • hda_ideal_VLE.py

  • hda_reaction.py

The state variables chosen for the property package are flows of component by phase, temperature and pressure. The components considered are: toluene, hydrogen, benzene and methane. Therefore, every stream has 8 flow variables, 1 temperature and 1 pressure variable.

1 Importing Modules#

1.1 Importing required Pyomo and IDAES components#

To construct a flowsheet, we will need several components from the Pyomo and IDAES package. Let us first import the following components from Pyomo:

  • Constraint (to write constraints)

  • Var (to declare variables)

  • ConcreteModel (to create the concrete model object)

  • Expression (to evaluate values as a function of variables defined in the model)

  • Objective (to define an objective function for optimization)

  • SolverFactory (to solve the problem)

  • TransformationFactory (to apply certain transformations)

  • Arc (to connect two unit models)

  • SequentialDecomposition (to initialize the flowsheet in a sequential mode)

For further details on these components, please refer to the Pyomo documentation: https://Pyomo.readthedocs.io/en/stable/

from pyomo.environ import (
    Constraint,
    Var,
    ConcreteModel,
    Expression,
    Objective,
    TransformationFactory,
    value,
)
from pyomo.network import Arc

From IDAES, we will be needing the FlowsheetBlock and the following unit models:

  • Feed

  • Mixer

  • Heater

  • StoichiometricReactor

  • Flash

  • Separator (splitter)

  • PressureChanger

  • Product

from idaes.core import FlowsheetBlock
from idaes.models.unit_models import (
    PressureChanger,
    Mixer,
    Separator as Splitter,
    Heater,
    StoichiometricReactor,
    Feed,
    Product,
)
Inline Exercise: Now, import the remaining unit models highlighted in blue above and run the cell using `Shift+Enter` after typing in the code.
# Todo: import flash model from idaes.models.unit_models
from idaes.models.unit_models import Flash

We will also be needing some utility tools to put together the flowsheet and calculate the degrees of freedom.

from idaes.models.unit_models.pressure_changer import ThermodynamicAssumption
from idaes.core.util.model_statistics import degrees_of_freedom

# Import idaes logger to set output levels
from idaes.core.solvers import get_solver

1.2 Importing required thermo and reaction package#

The final set of imports are to import the thermo and reaction package for the HDA process. We have created a custom thermo package that assumes Ideal Gas with support for VLE.

The reaction package here is very simple as we will be using only a StochiometricReactor and the reaction package consists of the stochiometric coefficients for the reaction and the parameter for the heat of reaction.

Let us import the following modules and they are in the same directory as this jupyter notebook:

  • hda_ideal_VLE as thermo_props
  • hda_reaction as reaction_props

from idaes.models.properties.modular_properties.base.generic_property import (
    GenericParameterBlock,
)
from idaes.models.properties.modular_properties.base.generic_reaction import (
    GenericReactionParameterBlock,
)
from idaes_examples.mod.hda.hda_ideal_VLE_modular import thermo_config
from idaes_examples.mod.hda.hda_reaction_modular import reaction_config

2 Constructing the Flowsheet#

We have now imported all the components, unit models, and property modules we need to construct a flowsheet. Let us create a ConcreteModel and add the flowsheet block.

m = ConcreteModel()
m.fs = FlowsheetBlock(dynamic=False)

We now need to add the property packages to the flowsheet. Unlike Module 1, where we only had a thermo property package, for this flowsheet we will also need to add a reaction property package.

m.fs.thermo_params = GenericParameterBlock(**thermo_config)
m.fs.reaction_params = GenericReactionParameterBlock(
    property_package=m.fs.thermo_params, **reaction_config
)

2.1 Adding Unit Models#

Let us start adding the unit models we have imported to the flowsheet. Here, we are adding the Feed (assigned a name I101 for Inlet), Mixer (assigned a name M101) and a Heater (assigned a name H101). Note that, all unit models need to be given a property package argument. In addition to that, there are several arguments depending on the unit model, please refer to the documentation for more details (https://idaes-pse.readthedocs.io/en/stable/reference_guides/model_libraries/generic/unit_models/index.html). For example, the Mixer unit model here must be specified the number of inlets that it will take in and the Heater can have specific settings enabled such as has_pressure_change or has_phase_equilibrium.

m.fs.I101 = Feed(property_package=m.fs.thermo_params)
m.fs.I102 = Feed(property_package=m.fs.thermo_params)

m.fs.M101 = Mixer(
    property_package=m.fs.thermo_params,
    num_inlets=3,
)

m.fs.H101 = Heater(
    property_package=m.fs.thermo_params,
    has_pressure_change=False,
    has_phase_equilibrium=True,
)
# Todo: Add reactor with the specifications above
m.fs.R101 = StoichiometricReactor(
    property_package=m.fs.thermo_params,
    reaction_package=m.fs.reaction_params,
    has_heat_of_reaction=True,
    has_heat_transfer=True,
    has_pressure_change=False,
)

Let us now add the Flash(assign the name F101) and pass the following arguments:

  • ”property_package”: m.fs.thermo_params
  • ”has_heat_transfer”: True
  • ”has_pressure_change”: False

m.fs.F101 = Flash(
    property_package=m.fs.thermo_params,
    has_heat_transfer=True,
    has_pressure_change=True,
)

Let us now add the Splitter(S101) with specific names for its output (purge and recycle), PressureChanger(C101) and the second Flash(F102).

m.fs.S101 = Splitter(
    property_package=m.fs.thermo_params,
    ideal_separation=False,
    outlet_list=["purge", "recycle"],
)


m.fs.C101 = PressureChanger(
    property_package=m.fs.thermo_params,
    compressor=True,
    thermodynamic_assumption=ThermodynamicAssumption.isothermal,
)

m.fs.F102 = Flash(
    property_package=m.fs.thermo_params,
    has_heat_transfer=True,
    has_pressure_change=True,
)

Last, we will add the three Product blocks (P101, P102, P103). We use Feed blocks and Product blocks for convenience with reporting stream summaries and consistency

m.fs.P101 = Product(property_package=m.fs.thermo_params)
m.fs.P102 = Product(property_package=m.fs.thermo_params)
m.fs.P103 = Product(property_package=m.fs.thermo_params)

2.2 Connecting Unit Models using Arcs#

We have now added all the unit models we need to the flowsheet. However, we have not yet specified how the units are to be connected. To do this, we will be using the Arc which is a Pyomo component that takes in two arguments: source and destination. Let us connect the outlet of the inlets (I101, I102) to the inlet of the mixer (M101) and outlet of the mixer to the inlet of the heater(H101).

m.fs.s01 = Arc(source=m.fs.I101.outlet, destination=m.fs.M101.inlet_1)
m.fs.s02 = Arc(source=m.fs.I102.outlet, destination=m.fs.M101.inlet_2)
m.fs.s03 = Arc(source=m.fs.M101.outlet, destination=m.fs.H101.inlet)
# Todo: Connect the H101 outlet to R101 inlet
m.fs.s04 = Arc(source=m.fs.H101.outlet, destination=m.fs.R101.inlet)

We will now be connecting the rest of the flowsheet as shown below. Notice how the outlet names are different for the flash tanks F101 and F102 as they have a vapor and a liquid outlet.

m.fs.s05 = Arc(source=m.fs.R101.outlet, destination=m.fs.F101.inlet)
m.fs.s06 = Arc(source=m.fs.F101.vap_outlet, destination=m.fs.S101.inlet)
m.fs.s07 = Arc(source=m.fs.F101.liq_outlet, destination=m.fs.F102.inlet)
m.fs.s08 = Arc(source=m.fs.S101.recycle, destination=m.fs.C101.inlet)
m.fs.s09 = Arc(source=m.fs.C101.outlet, destination=m.fs.M101.inlet_3)

Last we will connect the outlet streams to the inlets of the Product blocks (P101, P102, P103)

m.fs.s10 = Arc(source=m.fs.F102.vap_outlet, destination=m.fs.P101.inlet)
m.fs.s11 = Arc(source=m.fs.F102.liq_outlet, destination=m.fs.P102.inlet)
m.fs.s12 = Arc(source=m.fs.S101.purge, destination=m.fs.P103.inlet)

We have now connected the unit model block using the arcs. However, each of these arcs link to ports on the two unit models that are connected. In this case, the ports consist of the state variables that need to be linked between the unit models. Pyomo provides a convenient method to write these equality constraints for us between two ports and this is done as follows:

TransformationFactory("network.expand_arcs").apply_to(m)

2.3 Adding expressions to compute purity and operating costs#

In this section, we will add a few Expressions that allows us to evaluate the performance. Expressions provide a convenient way of calculating certain values that are a function of the variables defined in the model. For more details on Expressions, please refer to: https://pyomo.readthedocs.io/en/stable/explanation/modeling/network.html.

For this flowsheet, we are interested in computing the purity of the product Benzene stream (i.e. the mole fraction) and the operating cost which is a sum of the cooling and heating cost.

Let us first add an Expression to compute the mole fraction of benzene in the vap_outlet of F102 which is our product stream. Please note that the var flow_mol_phase_comp has the index - [time, phase, component]. As this is a steady-state flowsheet, the time index by default is 0. The valid phases are [“Liq”, “Vap”]. Similarly the valid component list is [“benzene”, “toluene”, “hydrogen”, “methane”].

m.fs.purity = Expression(
    expr=m.fs.F102.control_volume.properties_out[0].flow_mol_phase_comp[
        "Vap", "benzene"
    ]
    / (
        m.fs.F102.control_volume.properties_out[0].flow_mol_phase_comp["Vap", "benzene"]
        + m.fs.F102.control_volume.properties_out[0].flow_mol_phase_comp[
            "Vap", "toluene"
        ]
    )
)

Now, let us add an expression to compute the cooling cost assuming a cost of 0.212E-4 $/kW. Note that cooling utility is required for the reactor (R101) and the first flash (F101).

m.fs.cooling_cost = Expression(
    expr=0.212e-7 * (-m.fs.F101.heat_duty[0]) + 0.212e-7 * (-m.fs.R101.heat_duty[0])
)

Now, let us add an expression to compute the heating cost assuming the utility cost as follows:

  • 2.2E-4 dollars/kW for H101
  • 1.9E-4 dollars/kW for F102
Note that the heat duty is in units of watt (J/s).

m.fs.heating_cost = Expression(
    expr=2.2e-7 * m.fs.H101.heat_duty[0] + 1.9e-7 * m.fs.F102.heat_duty[0]
)

Let us now add an expression to compute the total operating cost per year which is basically the sum of the cooling and heating cost we defined above.

m.fs.operating_cost = Expression(
    expr=(3600 * 24 * 365 * (m.fs.heating_cost + m.fs.cooling_cost))
)

4 Specifying the Model#

4.1 Fixing feed conditions#

Let us first check how many degrees of freedom exist for this flowsheet using the degrees_of_freedom tool we imported earlier.

print(degrees_of_freedom(m))
29

We will now be fixing the toluene feed (I101) stream to the conditions shown in the flowsheet above. Please note that though this is a pure toluene feed, the remaining components are still assigned a very small non-zero value to help with convergence and initializing. We will be importing a function that will specify the inlet conditions for this example.

from idaes_examples.mod.hda.hda_flowsheet_extras import fix_inlet_states

tear_guesses = fix_inlet_states(m)

4.2 Fixing unit model specifications#

Now that we have fixed our inlet feed conditions, we will now be fixing the operating conditions for the unit models in the flowsheet. Let us set set the H101 outlet temperature to 600 K.

m.fs.H101.outlet.temperature.fix(600)

For the StoichiometricReactor, we have to define the conversion in terms of toluene. This requires us to create a new variable for specifying the conversion and adding a Constraint that defines the conversion with respect to toluene. The second degree of freedom for the reactor is to define the heat duty. In this case, let us assume the reactor to be adiabatic i.e. Q = 0.

m.fs.R101.conversion = Var(initialize=0.75, bounds=(0, 1))

m.fs.R101.conv_constraint = Constraint(
    expr=m.fs.R101.conversion
    * (m.fs.R101.control_volume.properties_in[0].flow_mol_phase_comp["Vap", "toluene"])
    == (
        m.fs.R101.control_volume.properties_in[0].flow_mol_phase_comp["Vap", "toluene"]
        - m.fs.R101.control_volume.properties_out[0].flow_mol_phase_comp[
            "Vap", "toluene"
        ]
    )
)

m.fs.R101.conversion.fix(0.75)
m.fs.R101.heat_duty.fix(0)

The Flash conditions for F101 can be set as follows.

m.fs.F101.vap_outlet.temperature.fix(325.0)
m.fs.F101.deltaP.fix(0)
m.fs.F102.vap_outlet.temperature.fix(375)
m.fs.F102.deltaP.fix(-200000)

Let us fix the purge split fraction to 20% and the outlet pressure of the compressor is set to 350000 Pa.

m.fs.S101.split_fraction[0, "purge"].fix(0.2)
m.fs.C101.outlet.pressure.fix(350000)
print(degrees_of_freedom(m))
0

5 Initializing the Model#

When a flowsheet contains a recycle loop, the outlet of a downstream unit becomes the inlet of an upstream unit, creating a cyclic dependency that prevents straightforward calculation of all stream conditions. The tear‐stream method is necessary because it “breaks” this loop: you select one recycle stream as the tear, assign it an initial guess, and then solve the rest of the flowsheet as if it were acyclic. Once the downstream units compute their outputs, you compare the calculated value of the torn stream to your initial guess and iteratively adjust until they coincide. Without tearing, the solver cannot establish a proper topological sequence or drive the recycle to convergence, making initialization—and ultimately steady‐state convergence—impossible.

It is important to determine the tear stream for a flowsheet which will be demonstrated below.

Currently, there are two methods of initializing a full flowsheet: using the sequential decomposition tool, or manually propagating through the flowsheet. The tear stream in this example will be the stream from the mixer to the heater since that is where the recycle stream first enters back into the main process.

First, we will highlight some helpful functions that are used in the initialization process.

def initialize_unit(unit):
    from idaes.core.util.exceptions import InitializationError
    import idaes.logger as idaeslog

    optarg = {
        "nlp_scaling_method": "user-scaling",
        "OF_ma57_automatic_scaling": "yes",
        "max_iter": 1000,
        "tol": 1e-8,
    }

    try:
        initializer = unit.default_initializer(solver_options=optarg)
        initializer.initialize(unit, output_level=idaeslog.INFO_LOW)
    except InitializationError:
        solver = get_solver(solver_options=optarg)
        solver.solve(unit)

This first function will take any unit model and can either initialize the model with its respective default initializer, or use a generic solver and the solve the current state of the unit model. Often times a direct initialization method will fail while a solving method will converge so having the option for both is helpful.

5.1 Sequential Decomposition#

This section will demonstrate how to use the built-in sequential decomposition tool to initialize our flowsheet. Sequential Decomposition is a tool from Pyomo where the documentation can be found here https://Pyomo.readthedocs.io/en/stable/explanation/modeling/network.html#sequential-decomposition

def automatic_propagation(m, tear_guesses):

    from pyomo.network import SequentialDecomposition

    seq = SequentialDecomposition()
    seq.options.select_tear_method = "heuristic"
    seq.options.tear_method = "Wegstein"
    seq.options.iterLim = 5

    # Using the SD tool
    G = seq.create_graph(m)
    heuristic_tear_set = seq.tear_set_arcs(G, method="heuristic")
    order = seq.calculation_order(G)

    # Pass the tear_guess to the SD tool
    seq.set_guesses_for(heuristic_tear_set[0].destination, tear_guesses)

    print(f"Tear Stream starts at: {heuristic_tear_set[0].destination.name}")

    for o in order:
        print(o[0].name)

    seq.run(m, initialize_unit)

We are now ready to initialize our flowsheet in a sequential mode. Note that we specifically set the iteration limit to be 5 as we are trying to use this tool only to get a good set of initial values such that IPOPT can then take over and solve this flowsheet for us. Uncomment this function call to run the automatic propagation method

# automatic_propagation(m, tear_guesses)

5.2 Manual Propagation Method#

This method uses a more direct approach to initialize the flowsheet, using the updated initializer method and propagating manually through the flowsheet and solving for the tear stream directly. Lets define the function that will help us manually propagate and step through the flowsheet

def manual_propagation(m, tear_guesses):
    from idaes.core.util.initialization import propagate_state

    print(f"The DOF is {degrees_of_freedom(m)} initially")
    m.fs.s03_expanded.deactivate()
    print(f"The DOF is {degrees_of_freedom(m)} after deactivating the tear stream")

    for k, v in tear_guesses.items():
        for k1, v1 in v.items():
            getattr(m.fs.s03.destination, k)[k1].fix(v1)

    print(f"The DOF is {degrees_of_freedom(m)} after setting the tear stream")

    optarg = {
        "nlp_scaling_method": "user-scaling",
        "OF_ma57_automatic_scaling": "yes",
        "max_iter": 300,
        # "tol": 1e-10,
    }

    solver = get_solver(solver_options=optarg)

    initialize_unit(m.fs.H101)  # Initialize Heater
    propagate_state(m.fs.s04)  # Establish connection between Heater and Reactor
    initialize_unit(m.fs.R101)  # Initialize Reactor
    propagate_state(
        m.fs.s05
    )  # Establish connection between Reactor and First Flash Unit
    initialize_unit(m.fs.F101)  # Initialize First Flash Unit
    propagate_state(
        m.fs.s06
    )  # Establish connection between First Flash Unit and Splitter
    propagate_state(
        m.fs.s07
    )  # Establish connection between First Flash Unit and Second Flash Unit
    initialize_unit(m.fs.S101)  # Initialize Splitter
    propagate_state(m.fs.s08)  # Establish connection between Splitter and Compressor
    initialize_unit(m.fs.C101)  # Initialize Compressor
    propagate_state(m.fs.s09)  # Establish connection between Compressor and Mixer
    initialize_unit(m.fs.I101)  # Initialize Toluene Inlet
    propagate_state(m.fs.s01)  # Establish connection between Toluene Inlet and Mixer
    initialize_unit(m.fs.I102)  # Initialize Hydrogen Inlet
    propagate_state(m.fs.s02)  # Establish connection between Hydrogen Inlet and Mixer
    initialize_unit(m.fs.M101)  # Initialize Mixer
    propagate_state(m.fs.s03)  # Establish connection between Mixer and Heater
    solver.solve(m.fs.F102)
    propagate_state(
        m.fs.s10
    )  # Establish connection between Second Flash Unit and Benzene Product
    propagate_state(
        m.fs.s11
    )  # Establish connection between Second Flash Unit and Toluene Product
    propagate_state(m.fs.s12)  # Establish connection between Splitter and Purge Product

    optarg = {
        "nlp_scaling_method": "user-scaling",
        "OF_ma57_automatic_scaling": "yes",
        "max_iter": 300,
        "tol": 1e-8,
    }
    solver = get_solver("ipopt_v2", options=optarg)
    solver.solve(m, tee=False)

    for k, v in tear_guesses.items():
        for k1, v1 in v.items():
            getattr(m.fs.H101.inlet, k)[k1].unfix()

    m.fs.s03_expanded.activate()
    print(
        f"The DOF is {degrees_of_freedom(m)} after unfixing the values and reactivating the tear stream"
    )

It will first show that the degrees of freedom is correctly at 0 before any streams are deactivated. Once the tear stream is deactivated though, the degrees of freedom will be 10. That means 10 variables will have to be defined with the tear guesses tear_guesses. Then each unit model can be initialized with our same helper function and then can propagate the corresponding connection to the following unit models. At the end, the whole flowsheet is solved, giving a much better chance for the recycle stream to be used correctly the flowsheet to converge.

manual_propagation(m, tear_guesses)
The DOF is 0 initially
The DOF is 10 after deactivating the tear stream
The DOF is 0 after setting the tear stream
The DOF is 0 after unfixing the values and reactivating the tear stream

6 Solving the Model#

We have now initialized the flowsheet. Lets set up some solving options before simulating the flowsheet. We want to specify the scaling method, number of iterations, and tolerance. More specific or advanced options can be found at the documentation for IPOPT https://coin-or.github.io/Ipopt/OPTIONS.html

optarg = {
    "nlp_scaling_method": "user-scaling",
    "OF_ma57_automatic_scaling": "yes",
    "max_iter": 1000,
    "tol": 1e-8,
}
# Create the solver object
solver = get_solver("ipopt_v2", options=optarg)

# Solve the model
results = solver.solve(m, tee=False)
Ipopt 3.13.2: linear_solver="ma57"
max_iter=1000
nlp_scaling_method="user-scaling"
tol=1e-08
option_file_name="C:\Users\Tanner\AppData\Local\Temp\tmpjzf4e172\unknown.31496.43076.opt"

Using option file "C:\Users\Tanner\AppData\Local\Temp\tmpjzf4e172\unknown.31496.43076.opt".


******************************************************************************
This program contains Ipopt, a library for large-scale nonlinear optimization.
 Ipopt is released as open source code under the Eclipse Public License (EPL).
         For more information visit http://projects.coin-or.org/Ipopt

This version of Ipopt was compiled from source code available at
    https://github.com/IDAES/Ipopt as part of the Institute for the Design of
    Advanced Energy Systems Process Systems Engineering Framework (IDAES PSE
    Framework) Copyright (c) 2018-2019. See https://github.com/IDAES/idaes-pse.

This version of Ipopt was compiled using HSL, a collection of Fortran codes
    for large-scale scientific computation.  All technical papers, sales and
    publicity material resulting from use of the HSL codes within IPOPT must
    contain the following acknowledgement:
        HSL, a collection of Fortran codes for large-scale scientific
        computation. See http://www.hsl.rl.ac.uk.
******************************************************************************

This is Ipopt version 3.13.2, running with linear solver ma57.

Number of nonzeros in equality constraint Jacobian...:      920
Number of nonzeros in inequality constraint Jacobian.:        0
Number of nonzeros in Lagrangian Hessian.............:      456

Total number of variables............................:      218
                     variables with only lower bounds:       56
                variables with lower and upper bounds:      155
                     variables with only upper bounds:        0
Total number of equality constraints.................:      218
Total number of inequality constraints...............:        0
        inequality constraints with only lower bounds:        0
   inequality constraints with lower and upper bounds:        0
        inequality constraints with only upper bounds:        0

iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
   0  0.0000000e+00 8.17e+03 1.00e+00  -1.0 0.00e+00    -  0.00e+00 0.00e+00   0
Reallocating memory for MA57: lfact (10193)
   1  0.0000000e+00 7.12e+03 2.01e+02  -1.0 4.35e+04    -  4.68e-02 1.56e-01h  1
   2  0.0000000e+00 7.12e+03 2.74e+02  -1.0 4.20e+05    -  9.82e-04 5.24e-04h  1
   3  0.0000000e+00 7.11e+03 2.06e+03  -1.0 4.26e+05    -  1.43e-05 1.49e-03f  1
   4  0.0000000e+00 7.09e+03 1.83e+03  -1.0 3.65e+04    -  5.78e-03 2.02e-03h  1
   5  0.0000000e+00 7.09e+03 1.82e+03  -1.0 1.61e+05    -  4.69e-04 1.23e-05h  1
   6r 0.0000000e+00 7.09e+03 9.99e+02   3.9 0.00e+00    -  0.00e+00 6.02e-08R  2
   7r 0.0000000e+00 6.68e+03 9.98e+02   3.9 7.11e+06    -  5.22e-04 1.63e-04f  1
   8r 0.0000000e+00 4.74e+04 9.98e+02   3.9 3.89e+06    -  2.84e-04 7.21e-04f  1
   9r 0.0000000e+00 4.77e+04 1.78e+04   3.9 1.73e+06    -  4.34e-02 5.81e-04f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
  10r 0.0000000e+00 4.55e+04 4.83e+04   3.2 4.71e+04    -  2.92e-01 4.47e-02f  1
  11r 0.0000000e+00 7.86e+04 3.48e+04   3.2 3.13e+02    -  4.51e-01 3.27e-01f  1
  12r 0.0000000e+00 4.85e+04 2.80e+04   3.2 5.27e+01    -  5.22e-01 2.61e-01f  1
  13r 0.0000000e+00 1.68e+05 1.76e+04   3.2 1.46e+02    -  8.30e-01 4.81e-01f  1
  14r 0.0000000e+00 8.02e+04 3.21e+03   3.2 1.42e+02    -  7.28e-01 1.00e+00f  1
  15r 0.0000000e+00 1.89e+05 6.22e+04   3.2 1.06e+02    -  4.83e-01 1.00e+00f  1
  16r 0.0000000e+00 1.77e+05 5.95e+04   3.2 1.08e+02   0.0 9.89e-02 6.28e-02h  1
  17r 0.0000000e+00 1.03e+05 3.98e+04   3.2 8.09e+01    -  7.18e-01 4.22e-01h  1
  18r 0.0000000e+00 8.88e+04 8.50e+04   3.2 1.59e+02    -  2.14e-01 1.36e-01h  1
  19r 0.0000000e+00 5.30e+04 1.71e+06   3.2 8.69e+01    -  3.02e-02 3.51e-01f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
  20r 0.0000000e+00 1.40e+04 2.20e+05   3.2 1.08e+01   2.2 6.72e-01 8.78e-01f  1
  21r 0.0000000e+00 9.61e+03 9.51e+04   2.5 1.06e+01   1.8 9.32e-01 4.27e-01f  1
  22r 0.0000000e+00 1.83e+04 3.37e+04   2.5 1.51e+01   1.3 1.00e+00 6.41e-01f  1
  23r 0.0000000e+00 3.76e+03 1.35e+03   2.5 4.90e+00   0.8 1.00e+00 1.00e+00f  1
  24r 0.0000000e+00 1.38e+05 2.72e+04   2.5 1.06e+02    -  3.87e-01 1.00e+00f  1
  25r 0.0000000e+00 1.06e+05 2.11e+04   2.5 5.11e+01   1.2 2.99e-01 2.35e-01h  1
  26r 0.0000000e+00 9.35e+04 1.83e+04   2.5 3.27e+01   2.6 6.94e-02 1.19e-01h  1
  27r 0.0000000e+00 9.23e+03 3.61e+03   2.5 2.06e+01   2.1 2.35e-01 1.00e+00h  1
  28r 0.0000000e+00 9.21e+03 6.60e+03   2.5 9.05e+02    -  1.37e-01 3.89e-03f  1
  29r 0.0000000e+00 6.69e+04 1.34e+04   2.5 1.58e+02    -  3.39e-01 2.01e-01f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
  30r 0.0000000e+00 7.42e+01 3.52e+04   2.5 3.26e-01   1.6 4.89e-01 1.00e+00f  1
  31r 0.0000000e+00 9.36e+03 1.66e+02   2.5 1.29e+01    -  1.00e+00 1.00e+00f  1
  32r 0.0000000e+00 9.30e+01 3.75e+00   2.5 1.85e+00    -  1.00e+00 1.00e+00h  1
  33r 0.0000000e+00 8.97e+03 9.16e+02   0.4 4.59e+01    -  7.14e-01 7.41e-01f  1
  34r 0.0000000e+00 7.19e+03 3.37e+02   0.4 2.34e+03    -  7.23e-01 6.68e-01f  1
  35r 0.0000000e+00 5.64e+03 5.42e+02   0.4 7.33e+02    -  6.81e-01 5.42e-01f  1
  36r 0.0000000e+00 1.10e+03 1.79e+02   0.4 2.38e-01   2.0 6.12e-01 8.67e-01f  1
  37r 0.0000000e+00 3.67e+02 1.23e+02   0.4 9.76e-02   2.5 9.90e-01 8.27e-01f  1
  38r 0.0000000e+00 2.93e+01 3.80e+02   0.4 2.55e-01   2.0 5.41e-01 1.00e+00f  1
  39r 0.0000000e+00 2.58e+01 7.52e+02   0.4 2.70e+00   2.4 4.08e-02 9.27e-02h  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
  40r 0.0000000e+00 1.03e+01 2.85e+02   0.4 8.74e-02   2.8 7.74e-01 6.78e-01f  1
  41r 0.0000000e+00 4.85e+00 8.16e+01   0.4 2.96e-01   2.4 7.04e-01 7.21e-01f  1
  42r 0.0000000e+00 1.15e+02 7.32e+01  -0.3 1.76e-01   1.9 7.14e-01 7.59e-01f  1
  43r 0.0000000e+00 2.10e+02 3.61e+02  -0.3 4.28e-01   1.4 7.77e-01 1.77e-01f  1
  44r 0.0000000e+00 1.51e+02 6.63e+01  -0.3 4.68e-01   0.9 1.00e+00 8.94e-01f  1
  45r 0.0000000e+00 9.19e+00 1.85e+01  -0.3 1.49e-01   1.3 1.00e+00 1.00e+00f  1
  46r 0.0000000e+00 5.62e+00 5.71e+01  -0.3 1.05e+00   0.9 4.17e-01 3.89e-01h  1
  47r 0.0000000e+00 1.98e+01 3.77e+01  -0.3 7.16e-01   0.4 1.00e+00 1.00e+00f  1
  48r 0.0000000e+00 1.82e+01 6.27e+02  -0.3 9.97e+00  -0.1 4.98e-02 8.93e-02h  1
  49r 0.0000000e+00 1.72e+01 5.84e+02  -0.3 1.90e+00   1.2 4.67e-02 4.96e-02f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
  50r 0.0000000e+00 1.16e+01 7.04e+02  -0.3 1.10e+00   0.8 6.53e-01 4.12e-01f  1
  51r 0.0000000e+00 1.08e+01 6.28e+02  -0.3 2.94e-01   2.1 3.66e-01 6.78e-02f  1
  52r 0.0000000e+00 6.69e+00 7.88e+02  -0.3 4.74e-01   1.6 5.14e-01 6.14e-01f  1
  53r 0.0000000e+00 6.34e+00 1.39e+03  -0.3 6.49e-02   2.9 1.00e+00 5.17e-02h  1
  54r 0.0000000e+00 5.30e-01 7.34e+02  -0.3 1.91e-01   2.5 5.47e-01 1.00e+00f  1
  55r 0.0000000e+00 5.30e-01 8.37e+02  -0.3 1.09e-01   2.9 1.00e+00 1.86e-01f  1
  56r 0.0000000e+00 6.24e-01 3.25e+02  -0.3 2.18e-01   2.4 7.22e-01 7.62e-01f  1
  57r 0.0000000e+00 5.84e-01 1.25e+03  -0.3 1.10e-01   2.8 3.44e-01 6.58e-02f  1
  58r 0.0000000e+00 1.84e+00 9.86e+02  -0.3 6.56e-01   2.4 4.76e-01 1.66e-01f  1
  59r 0.0000000e+00 9.98e+00 7.31e+02  -0.3 2.98e-01   1.9 3.81e-01 1.96e-01f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
  60r 0.0000000e+00 2.33e+01 6.73e+02  -0.3 2.93e+00   1.4 1.64e-03 2.62e-02f  1
  61r 0.0000000e+00 2.09e+02 6.21e+02  -0.3 3.68e+00   1.8 2.80e-02 4.31e-02f  1
  62r 0.0000000e+00 2.06e+02 4.81e+02  -0.3 4.45e-01   1.4 4.03e-01 1.68e-02h  1
  63r 0.0000000e+00 1.65e+02 3.87e+02  -0.3 1.11e+00   0.9 1.92e-02 1.96e-01f  1
  64r 0.0000000e+00 1.50e+02 3.08e+02  -0.3 1.85e+00   0.4 2.04e-01 2.07e-01f  1
  65r 0.0000000e+00 1.42e+02 2.98e+02  -0.3 1.20e+00   1.7 3.04e-02 5.07e-02f  1
  66r 0.0000000e+00 7.46e+01 5.29e+02  -0.3 5.27e-02   2.2 9.16e-01 4.78e-01f  1
  67r 0.0000000e+00 5.57e+01 4.82e+02  -0.3 1.72e-01   1.7 9.61e-01 2.54e-01f  1
  68r 0.0000000e+00 2.98e+00 1.46e+02  -0.3 2.22e-01   1.2 1.00e+00 9.79e-01f  1
  69r 0.0000000e+00 6.01e+00 2.34e+02  -0.3 4.07e-01   0.7 5.24e-01 4.28e-01f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
  70r 0.0000000e+00 1.07e+02 5.44e+02  -0.3 1.11e+00   0.3 7.04e-01 1.00e+00f  1
  71r 0.0000000e+00 2.25e+01 7.39e+02  -0.3 1.39e-01   1.6 3.98e-01 7.90e-01f  1
  72r 0.0000000e+00 2.07e+01 7.10e+02  -0.3 6.94e-01   1.1 2.42e-01 8.16e-02f  1
  73r 0.0000000e+00 4.67e+01 3.79e+02  -0.3 3.86e-01   0.6 6.91e-01 1.00e+00f  1
  74r 0.0000000e+00 7.77e+01 2.65e+02  -0.3 9.27e-01   0.2 4.94e-01 4.38e-01h  1
  75r 0.0000000e+00 7.62e+01 2.61e+02  -0.3 3.81e+00   0.6 2.20e-02 3.06e-02f  2
  76r 0.0000000e+00 5.81e+01 2.16e+02  -0.3 1.57e-01   1.9 3.64e-01 2.38e-01f  1
  77r 0.0000000e+00 1.97e+01 1.33e+02  -0.3 1.61e-01   1.4 1.00e+00 6.92e-01f  1
  78r 0.0000000e+00 4.47e+01 2.11e+02  -0.3 8.61e-01   1.0 7.64e-01 4.57e-01f  1
  79r 0.0000000e+00 2.75e+02 1.68e+02  -0.3 2.89e+00   0.5 1.90e-01 1.97e-01f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
  80r 0.0000000e+00 1.57e+01 5.93e+01  -0.3 1.06e-01   1.8 1.00e+00 1.00e+00f  1
  81r 0.0000000e+00 4.52e+00 7.54e+00  -0.3 3.54e-01   1.3 1.00e+00 1.00e+00h  1
  82r 0.0000000e+00 1.51e+01 6.83e+01  -1.0 1.48e-01   1.8 6.09e-01 6.86e-01f  1
  83r 0.0000000e+00 9.47e+00 2.62e+01  -1.0 5.68e-02   2.2 1.00e+00 8.94e-01f  1
  84r 0.0000000e+00 2.30e+01 2.17e+01  -1.0 4.35e-01   1.7 1.62e-01 1.76e-01f  1
  85r 0.0000000e+00 1.67e+01 1.83e+02  -1.0 8.29e-02   2.1 1.00e+00 4.80e-01f  1
  86r 0.0000000e+00 5.83e+01 5.12e+02  -1.0 3.20e-01   1.7 6.90e-01 3.58e-01f  1
  87r 0.0000000e+00 1.07e+02 1.19e+03  -1.0 1.30e+01   1.2 2.99e-02 1.02e-02f  1
  88r 0.0000000e+00 1.05e+02 3.60e+02  -1.0 2.66e+00   0.7 1.00e+00 2.59e-02f  1
  89r 0.0000000e+00 4.13e+02 4.81e+02  -1.0 7.27e+00   0.2 1.00e+00 4.22e-01f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
  90r 0.0000000e+00 4.12e+02 4.73e+02  -1.0 7.91e+01  -0.3 2.91e-03 3.88e-03f  1
  91r 0.0000000e+00 3.90e+02 5.09e+02  -1.0 1.05e+00   1.1 1.00e+00 5.42e-02f  1
  92r 0.0000000e+00 6.42e+02 1.96e+01  -1.0 3.03e+00   0.6 1.00e+00 1.00e+00f  1
  93r 0.0000000e+00 1.54e+03 5.81e+02  -1.0 1.25e+01   0.1 9.40e-02 6.11e-01f  1
  94r 0.0000000e+00 1.52e+03 5.61e+02  -1.0 9.43e+00   0.5 1.64e-02 8.32e-03f  1
  95r 0.0000000e+00 1.38e+03 3.41e+02  -1.0 2.32e+01   0.1 3.92e-01 1.27e-01f  1
  96r 0.0000000e+00 1.38e+03 3.72e+02  -1.0 5.95e+00   0.5 3.61e-01 4.34e-03h  1
  97r 0.0000000e+00 1.28e+03 2.01e+02  -1.0 1.09e+01   0.0 4.23e-02 1.05e-01f  1
  98r 0.0000000e+00 9.89e+02 4.20e+02  -1.0 4.94e+00   0.4 3.76e-02 2.66e-01f  1
  99r 0.0000000e+00 8.96e+02 4.65e+02  -1.0 2.96e+00   0.9 2.83e-01 9.51e-02f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
 100r 0.0000000e+00 8.38e+02 1.08e+03  -1.0 3.79e+00   0.4 5.07e-01 6.96e-02h  1
 101r 0.0000000e+00 1.19e+03 6.44e+02  -1.0 9.81e+00  -0.1 3.71e-01 3.57e-01f  1
 102r 0.0000000e+00 1.19e+03 5.28e+02  -1.0 1.78e+01  -0.6 7.44e-03 6.07e-02f  1
 103r 0.0000000e+00 1.02e+03 3.34e+02  -1.0 2.10e+01  -1.0 4.17e-02 1.21e-01f  1
 104r 0.0000000e+00 9.98e+02 8.36e+02  -1.0 8.35e+00  -0.6 4.71e-01 1.95e-02h  1
 105r 0.0000000e+00 6.83e+02 6.45e+02  -1.0 1.30e+01  -1.1 5.71e-02 2.84e-01f  1
 106r 0.0000000e+00 8.65e+02 1.14e+03  -1.0 8.25e+00  -0.7 1.61e-01 6.21e-01f  1
 107r 0.0000000e+00 7.16e+02 9.70e+02  -1.0 4.28e+00  -0.2 1.93e-01 1.57e-01h  1
 108r 0.0000000e+00 6.31e+02 8.28e+02  -1.0 6.49e+00  -0.7 6.61e-02 1.81e-01h  1
 109r 0.0000000e+00 4.55e+02 1.19e+03  -1.0 2.95e+00  -0.3 1.88e-01 8.74e-01h  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
 110r 0.0000000e+00 1.83e+02 4.89e+02  -1.0 1.15e-01   1.9 6.47e-01 5.98e-01h  1
 111r 0.0000000e+00 3.73e+01 1.31e+02  -1.0 8.98e-02   1.5 1.00e+00 7.95e-01h  1
 112r 0.0000000e+00 7.86e+00 8.49e+01  -1.0 7.76e-02   1.9 9.93e-01 7.86e-01h  1
 113r 0.0000000e+00 2.86e+00 1.42e+01  -1.0 1.00e-01   1.4 1.00e+00 1.00e+00h  1
 114r 0.0000000e+00 2.89e+00 2.59e+00  -1.0 3.00e-01   0.9 1.00e+00 1.00e+00f  1
 115r 0.0000000e+00 5.21e+01 6.10e+00  -1.0 8.19e-01   0.5 1.00e+00 1.00e+00f  1
 116r 0.0000000e+00 1.36e+02 4.46e+02  -1.0 2.16e+00  -0.0 2.00e-01 5.43e-01f  1
 117r 0.0000000e+00 2.35e+01 1.35e+02  -1.0 3.81e-02   2.2 6.64e-01 8.26e-01h  1
 118r 0.0000000e+00 1.45e+01 6.04e+02  -1.0 4.14e-01   1.7 1.30e-01 6.39e-01f  1
 119r 0.0000000e+00 7.09e+01 1.51e+03  -1.0 1.49e+00   2.2 9.44e-02 3.42e-01f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
 120r 0.0000000e+00 6.91e+01 1.74e+03  -1.0 3.65e-01   2.6 2.02e-01 2.71e-02f  1
 121r 0.0000000e+00 6.89e+01 1.47e+03  -1.0 1.29e+00   2.1 8.23e-04 3.05e-03f  1
 122r 0.0000000e+00 6.81e+01 2.26e+03  -1.0 1.01e+00   1.6 5.80e-02 1.11e-02f  1
 123r 0.0000000e+00 5.53e+01 8.13e+02  -1.0 9.48e-01   1.2 1.55e-02 1.67e-01f  1
 124r 0.0000000e+00 4.39e+01 6.51e+02  -1.0 6.05e-01   0.7 2.20e-01 1.95e-01f  1
 125r 0.0000000e+00 3.95e+01 5.85e+02  -1.0 1.15e+00   0.2 4.00e-01 9.82e-02f  1
 126r 0.0000000e+00 5.31e+01 2.97e+02  -1.0 1.75e+00  -0.3 4.11e-03 3.31e-01f  1
 127r 0.0000000e+00 5.22e+01 2.92e+02  -1.0 1.49e+00   0.2 4.34e-02 1.79e-02h  1
 128r 0.0000000e+00 5.15e+01 2.70e+02  -1.0 9.28e+00  -0.3 1.26e-01 3.54e-02f  1
 129r 0.0000000e+00 2.87e+01 8.41e+01  -1.0 1.20e+00   0.1 9.87e-01 7.24e-01h  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
 130r 0.0000000e+00 5.74e+01 7.12e+01  -1.0 2.02e+00  -0.4 1.00e+00 3.08e-01h  1
 131r 0.0000000e+00 1.40e+03 2.57e+01  -1.0 6.81e+00  -0.9 6.75e-01 8.23e-01f  1
 132r 0.0000000e+00 1.35e+03 2.54e+01  -1.0 1.62e+01  -1.3 1.03e-01 3.34e-02h  1
 133r 0.0000000e+00 1.23e+03 3.95e+01  -1.0 7.83e+00  -0.9 7.40e-01 6.60e-01h  1
 134r 0.0000000e+00 1.19e+03 3.30e+02  -1.0 1.81e+01  -1.4 3.93e-01 5.73e-02h  1
 135r 0.0000000e+00 1.18e+03 3.16e+02  -1.0 1.12e+02  -1.9 2.06e-02 2.35e-02h  1
 136r 0.0000000e+00 3.04e+03 1.74e+02  -1.0 1.66e+01  -1.4 6.20e-01 5.67e-01h  1
 137r 0.0000000e+00 3.10e+03 7.86e+01  -1.0 4.18e+02  -1.9 4.32e-04 5.13e-03f  1
 138r 0.0000000e+00 3.15e+03 8.59e+01  -1.0 2.08e+01  -1.5 6.10e-01 3.20e-01f  1
 139r 0.0000000e+00 6.79e+03 9.92e+01  -1.0 6.29e+01  -2.0 5.13e-01 5.43e-01h  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
 140r 0.0000000e+00 1.99e+03 2.93e+02  -1.0 1.82e-01   1.2 9.52e-01 7.07e-01h  1
 141r 0.0000000e+00 4.14e-01 6.19e+01  -1.0 4.03e-02   1.6 9.05e-01 1.00e+00h  1
 142r 0.0000000e+00 5.40e-01 1.00e+01  -1.0 7.66e-02   1.1 1.00e+00 1.00e+00h  1
 143r 0.0000000e+00 1.87e+00 7.98e-01  -1.0 1.80e-01   0.6 1.00e+00 1.00e+00h  1
 144r 0.0000000e+00 1.04e+01 9.40e-01  -1.0 4.36e-01   0.2 1.00e+00 1.00e+00h  1
 145r 0.0000000e+00 3.19e+01 4.06e+00  -1.0 9.89e-01  -0.3 1.00e+00 1.00e+00h  1
 146r 0.0000000e+00 2.28e+02 6.27e+00  -1.0 3.02e+00  -0.8 1.00e+00 9.35e-01h  1
 147r 0.0000000e+00 7.58e+02 1.93e+02  -1.0 1.17e+01  -1.3 1.00e+00 5.43e-01h  1
 148r 0.0000000e+00 7.62e+02 6.49e+02  -1.0 5.17e+01  -1.7 4.76e-01 9.10e-02h  1
 149r 0.0000000e+00 5.18e+03 4.98e+02  -1.0 3.35e+02  -2.2 2.78e-02 1.43e-01f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
 150r 0.0000000e+00 5.19e+03 7.41e+02  -1.0 1.23e+03  -2.7 1.42e-02 2.95e-02f  1
 151r 0.0000000e+00 6.87e+03 9.57e+02  -1.0 1.12e+03  -2.3 1.35e-03 3.17e-02f  1
 152r 0.0000000e+00 6.85e+03 8.50e+02  -1.0 1.37e+02  -1.8 3.99e-01 4.06e-03h  1
 153r 0.0000000e+00 6.72e+03 9.57e+02  -1.0 4.39e+02  -2.3 1.18e-01 1.81e-02h  1
 154r 0.0000000e+00 6.40e+03 1.13e+03  -1.0 1.07e+04    -  1.49e-01 4.89e-02f  1
 155r 0.0000000e+00 5.62e+03 9.92e+02  -1.0 9.94e+03    -  1.20e-01 1.22e-01f  1
 156r 0.0000000e+00 5.58e+03 9.64e+02  -1.0 3.35e+03    -  1.73e-03 7.88e-03h  1
 157r 0.0000000e+00 5.58e+03 9.64e+02  -1.0 2.55e+05  -2.8 2.32e-05 2.33e-05f  1
 158r 0.0000000e+00 5.55e+03 1.02e+03  -1.0 4.99e+02  -3.3 2.90e-02 5.00e-03h  1
 159r 0.0000000e+00 4.91e+03 9.36e+02  -1.0 1.06e+02  -3.8 1.63e-01 1.16e-01h  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
 160r 0.0000000e+00 4.19e+03 7.63e+02  -1.0 1.75e+02  -2.4 1.15e-02 1.48e-01h  1
 161r 0.0000000e+00 4.18e+03 7.94e+02  -1.0 4.34e+02  -2.9 3.97e-02 1.46e-03h  1
 162r 0.0000000e+00 2.31e+03 3.53e+02  -1.0 5.33e+03    -  1.01e-01 4.63e-01h  1
 163r 0.0000000e+00 2.04e+03 5.97e+02  -1.0 3.18e+03    -  1.73e-01 1.20e-01h  1
 164r 0.0000000e+00 1.26e+03 3.46e+02  -1.0 3.03e+03    -  4.18e-01 4.11e-01h  1
 165r 0.0000000e+00 1.24e+03 3.20e+02  -1.0 2.32e+03    -  1.97e-03 1.44e-02h  1
 166r 0.0000000e+00 1.17e+03 4.22e+02  -1.0 2.31e+03    -  3.37e-01 5.42e-02h  1
 167r 0.0000000e+00 1.04e+03 6.87e+02  -1.0 1.74e+03    -  9.08e-01 1.16e-01h  1
 168r 0.0000000e+00 3.09e+02 4.10e+02  -1.0 1.64e+03    -  4.91e-01 8.06e-01h  1
 169r 0.0000000e+00 3.95e+02 5.59e+02  -1.0 2.46e+03    -  1.49e-02 2.83e-02f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
 170r 0.0000000e+00 2.99e+02 7.69e+01  -1.0 4.31e+02    -  8.88e-01 1.00e+00f  1
 171r 0.0000000e+00 1.65e+02 3.84e+02  -1.0 4.91e-01   0.2 1.00e+00 4.50e-01h  1
 172r 0.0000000e+00 1.87e-01 5.61e+01  -1.0 2.62e+00    -  1.00e+00 1.00e+00h  1
 173r 0.0000000e+00 1.87e-01 2.09e+01  -1.0 8.62e-01    -  1.00e+00 1.00e+00h  1
 174r 0.0000000e+00 1.87e-01 1.77e+00  -1.0 9.94e-02    -  1.00e+00 1.00e+00h  1
 175r 0.0000000e+00 7.96e+00 5.39e+00  -1.7 1.21e+01    -  1.00e+00 9.98e-01f  1
 176r 0.0000000e+00 2.31e+02 1.33e+02  -1.7 1.26e+04    -  3.42e-01 1.24e-01f  1
 177r 0.0000000e+00 2.38e+02 3.38e+02  -1.7 1.05e+04    -  6.40e-01 3.94e-02f  1
 178r 0.0000000e+00 2.37e+02 3.35e+02  -1.7 7.28e+03    -  9.97e-02 1.97e-03h  1
 179r 0.0000000e+00 2.11e+02 3.90e+02  -1.7 3.85e+03    -  5.39e-01 1.07e-01h  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
 180r 0.0000000e+00 3.16e+02 4.17e+02  -1.7 3.82e+03    -  1.28e-02 9.51e-01f  1
 181r 0.0000000e+00 1.07e+03 2.19e+02  -1.7 3.96e+03    -  4.42e-01 4.76e-01h  1
 182r 0.0000000e+00 1.05e+03 2.05e+02  -1.7 1.41e+03    -  4.64e-01 2.04e-02h  1
 183r 0.0000000e+00 4.56e+02 1.83e+02  -1.7 1.51e+00  -0.2 4.60e-01 5.66e-01h  1
 184r 0.0000000e+00 4.46e+02 4.37e+02  -1.7 1.26e+03    -  5.85e-01 2.25e-02h  1
 185r 0.0000000e+00 3.62e+02 5.95e+02  -1.7 1.17e+03    -  3.15e-01 2.25e-01h  1
 186r 0.0000000e+00 3.47e+02 2.11e+02  -1.7 6.54e-01  -0.7 4.90e-01 4.12e-02h  1
 187r 0.0000000e+00 2.68e+02 4.10e+02  -1.7 9.73e+02    -  1.74e-01 2.84e-01h  1
 188r 0.0000000e+00 2.32e+02 7.23e+02  -1.7 1.37e+03    -  7.54e-01 3.83e-01h  1
 189r 0.0000000e+00 2.25e+02 1.95e+02  -1.7 7.05e+02    -  6.62e-01 3.01e-02h  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
 190r 0.0000000e+00 1.61e+02 3.01e+02  -1.7 4.09e+02    -  1.00e+00 3.03e-01h  1
 191r 0.0000000e+00 9.61e-01 6.13e-01  -1.7 6.25e+01    -  1.00e+00 1.00e+00h  1
 192r 0.0000000e+00 5.56e-02 3.83e-02  -1.7 6.74e+00    -  1.00e+00 1.00e+00h  1
 193r 0.0000000e+00 3.80e-01 2.41e+01  -3.9 2.90e+01    -  9.11e-01 7.83e-01f  1
 194r 0.0000000e+00 2.60e+00 2.67e+01  -3.9 4.29e+04    -  7.72e-02 5.58e-02f  1
 195r 0.0000000e+00 2.60e+00 3.11e+02  -3.9 2.37e+04    -  6.85e-02 1.39e-05h  1
 196r 0.0000000e+00 2.91e+00 2.90e+02  -3.9 3.11e+03    -  1.53e-01 3.55e-02f  1
 197r 0.0000000e+00 1.00e+01 2.50e+02  -3.9 2.55e+03    -  1.71e-01 1.67e-01f  1
 198r 0.0000000e+00 2.87e+01 3.43e+02  -3.9 2.12e+03    -  7.25e-02 3.31e-01f  1
 199r 0.0000000e+00 4.82e+01 9.78e+02  -3.9 1.55e+03    -  5.38e-02 5.45e-01f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
 200r 0.0000000e+00 1.40e+02 6.78e+02  -3.9 3.21e+04    -  2.94e-02 1.88e-01f  1
 201r 0.0000000e+00 1.38e+02 6.81e+02  -3.9 2.53e+04    -  4.87e-02 1.94e-02h  1
 202r 0.0000000e+00 1.43e+02 6.07e+02  -3.9 2.53e+04    -  2.36e-02 7.13e-02f  1
 203r 0.0000000e+00 2.17e+02 4.73e+02  -3.9 2.36e+04    -  9.43e-02 1.98e-01f  1
 204r 0.0000000e+00 2.11e+02 5.82e+02  -3.9 1.89e+04    -  5.56e-01 6.52e-02h  1
 205r 0.0000000e+00 2.11e+02 7.24e+02  -3.9 1.33e+04    -  8.62e-01 3.50e-03h  1
 206r 0.0000000e+00 1.88e+02 8.64e+02  -3.9 3.55e+02    -  1.00e+00 1.09e-01h  1
 207r 0.0000000e+00 3.04e+01 5.52e+02  -3.9 3.16e+02    -  2.11e-01 8.62e-01h  1
 208r 0.0000000e+00 7.41e-01 1.91e+00  -3.9 4.38e+01    -  1.00e+00 9.97e-01h  1
 209r 0.0000000e+00 7.13e-01 6.42e+02  -3.9 1.50e-01    -  7.36e-02 3.78e-02h  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
 210r 0.0000000e+00 5.37e-03 4.03e+02  -3.9 1.44e-01    -  5.86e-01 1.00e+00h  1
 211r 0.0000000e+00 5.37e-03 2.68e-04  -3.9 2.52e-03    -  1.00e+00 1.00e+00h  1
 212r 0.0000000e+00 5.37e-03 5.94e+01  -5.9 3.33e-01    -  1.00e+00 6.67e-01f  1
 213r 0.0000000e+00 4.67e+01 1.51e+02  -5.9 2.60e+04    -  1.03e-01 2.28e-02f  1
 214r 0.0000000e+00 4.67e+01 8.78e+02  -5.9 4.07e+02    -  8.75e-01 9.73e-07h  2
 215r 0.0000000e+00 1.68e+01 2.98e+02  -5.9 1.59e-01    -  8.59e-01 6.40e-01h  1
 216r 0.0000000e+00 2.59e+00 5.90e+01  -5.9 3.25e-02    -  1.00e+00 8.46e-01h  1
 217r 0.0000000e+00 1.57e-04 3.42e-04  -5.9 3.65e-03    -  1.00e+00 1.00e+00h  1
 218r 0.0000000e+00 2.69e-09 6.22e-06  -5.9 2.22e-04    -  1.00e+00 1.00e+00h  1

Number of Iterations....: 218

                                   (scaled)                 (unscaled)
Objective...............:   0.0000000000000000e+00    0.0000000000000000e+00
Dual infeasibility......:   0.0000000000000000e+00    0.0000000000000000e+00
Constraint violation....:   2.6943874215090702e-09    2.6943874215090702e-09
Complementarity.........:   0.0000000000000000e+00    0.0000000000000000e+00
Overall NLP error.......:   2.6943874215090702e-09    2.6943874215090702e-09


Number of objective function evaluations             = 225
Number of objective gradient evaluations             = 8
Number of equality constraint evaluations            = 225
Number of inequality constraint evaluations          = 0
Number of equality constraint Jacobian evaluations   = 220
Number of inequality constraint Jacobian evaluations = 0
Number of Lagrangian Hessian evaluations             = 218
Total CPU secs in IPOPT (w/o function evaluations)   =      0.347
Total CPU secs in NLP function evaluations           =      0.036

EXIT: Optimal Solution Found.

7 Analyze the results#

If the IDAES UI package was installed with the idaes-pse installation or installed separately, you can run the flowsheet visualizer to see a full diagram of the full process that is generated and displayed on a browser window.

Otherwise, we can run the m.fs.report() method to see a full summary of the solved flowsheet. It is recommended to adjust the width of the output as much as possible for the cleanest display.

m.fs.report()
====================================================================================
Flowsheet : fs                                                             Time: 0.0
------------------------------------------------------------------------------------
    Stream Table
                                                Units          s01        s02        s03        s04        s05        s06        s07        s08        s09        s10        s11        s12   
    Total Molar Flowrate Liq                 mole / second    0.30001 2.0000e-05    0.34190 1.6073e-09 5.7340e-09 1.0000e-08    0.26712 1.1139e-06 1.1143e-06 1.0000e-08   0.094878 2.7856e-07
    Total Molar Flowrate Vap                 mole / second 4.0000e-05    0.32002     1.6901     2.0320     2.0320     1.7648 1.0000e-08     1.4119     1.4119    0.17224 1.0000e-08    0.35297
    Total Mole Fraction ('Liq', 'benzene')   dimensionless 3.3332e-05    0.50000    0.22733    0.13374    0.63390    0.76595    0.76595    0.76595    0.76595    0.66001    0.66001    0.76595
    Total Mole Fraction ('Liq', 'toluene')   dimensionless    0.99997    0.50000    0.77267    0.86626    0.36610    0.23405    0.23405    0.23405    0.23405    0.33999    0.33999    0.23405
    Total Mole Fraction ('Vap', 'benzene')   dimensionless    0.25000 3.1248e-05   0.024624   0.058732    0.17408   0.084499   0.084499   0.084499   0.084499    0.82430    0.82430   0.084499
    Total Mole Fraction ('Vap', 'toluene')   dimensionless    0.25000 3.1248e-05   0.028601    0.15380   0.038450  0.0088437  0.0088437  0.0088435  0.0088435    0.17570    0.17570  0.0088435
    Total Mole Fraction ('Vap', 'hydrogen')  dimensionless    0.25000    0.93744    0.33283    0.27683    0.16148    0.18592    0.18592    0.18592    0.18592 1.7561e-08 1.7561e-08    0.18592
    Total Mole Fraction ('Vap', 'methane')   dimensionless    0.25000   0.062496    0.61394    0.51064    0.62599    0.72074    0.72074    0.72074    0.72074 4.3265e-08 4.3265e-08    0.72074
    Temperature                                     kelvin     303.20     303.20     324.51     600.00     771.86     325.00     325.00     325.00     325.00     375.00     375.00     325.00
    Pressure                                        pascal 3.5000e+05 3.5000e+05 3.5000e+05 3.5000e+05 3.5000e+05 3.5000e+05 3.5000e+05 3.5000e+05 3.5000e+05 1.5000e+05 1.5000e+05 3.5000e+05
====================================================================================

What is the total operating cost?

print("operating cost = $", value(m.fs.operating_cost))
operating cost = $ 419008.281895999

For this operating cost, what is the amount of benzene we are able to produce and what purity we are able to achieve? We can look at a specific unit models stream table with the same report() method.

m.fs.F102.report()

print()
print("benzene purity = ", value(m.fs.purity))
====================================================================================
Unit : fs.F102                                                             Time: 0.0
------------------------------------------------------------------------------------
    Unit Performance

    Variables: 

    Key             : Value       : Units  : Fixed : Bounds
          Heat Duty :      7346.7 :   watt : False : (None, None)
    Pressure Change : -2.0000e+05 : pascal :  True : (None, None)

------------------------------------------------------------------------------------
    Stream Table
                                                 Units         Inlet    Vapor Outlet  Liquid Outlet
    Total Molar Flowrate Liq                  mole / second    0.26712            -             -  
    Total Molar Flowrate Vap                  mole / second 1.0000e-08            -             -  
    Total Mole Fraction ('Liq', 'benzene')    dimensionless    0.76595            -             -  
    Total Mole Fraction ('Liq', 'toluene')    dimensionless    0.23405            -             -  
    Total Mole Fraction ('Vap', 'benzene')    dimensionless   0.084499            -             -  
    Total Mole Fraction ('Vap', 'toluene')    dimensionless  0.0088437            -             -  
    Total Mole Fraction ('Vap', 'hydrogen')   dimensionless    0.18592            -             -  
    Total Mole Fraction ('Vap', 'methane')    dimensionless    0.72074            -             -  
    Temperature                                      kelvin     325.00            -             -  
    Pressure                                         pascal 3.5000e+05            -             -  
    flow_mol_phase Liq                        mole / second          -   1.0000e-08      0.094878  
    flow_mol_phase Vap                        mole / second          -      0.17224    1.0000e-08  
    mole_frac_phase_comp ('Liq', 'benzene')   dimensionless          -      0.66001       0.66001  
    mole_frac_phase_comp ('Liq', 'toluene')   dimensionless          -      0.33999       0.33999  
    mole_frac_phase_comp ('Vap', 'benzene')   dimensionless          -      0.82430       0.82430  
    mole_frac_phase_comp ('Vap', 'toluene')   dimensionless          -      0.17570       0.17570  
    mole_frac_phase_comp ('Vap', 'hydrogen')  dimensionless          -   1.7561e-08    1.7561e-08  
    mole_frac_phase_comp ('Vap', 'methane')   dimensionless          -   4.3265e-08    4.3265e-08  
    temperature                                      kelvin          -       375.00        375.00  
    pressure                                         pascal          -   1.5000e+05    1.5000e+05  
====================================================================================

benzene purity =  0.8242963450787166

Next, let’s look at how much benzene we are losing with the light gases out of F101. IDAES has tools for creating stream tables based on the Arcs and/or Ports in a flowsheet. Let us create and print a simple stream table showing the stream leaving the reactor and the vapor stream from F101.

from idaes.core.util.tables import (
    create_stream_table_dataframe,
    stream_table_dataframe_to_string,
)

st = create_stream_table_dataframe({"Reactor": m.fs.s05, "Light Gases": m.fs.s06})
print(stream_table_dataframe_to_string(st))
                                            Units        Reactor   Light Gases
Total Molar Flowrate Liq                 mole / second 5.7340e-09  1.0000e-08 
Total Molar Flowrate Vap                 mole / second     2.0320      1.7648 
Total Mole Fraction ('Liq', 'benzene')   dimensionless    0.63390     0.76595 
Total Mole Fraction ('Liq', 'toluene')   dimensionless    0.36610     0.23405 
Total Mole Fraction ('Vap', 'benzene')   dimensionless    0.17408    0.084499 
Total Mole Fraction ('Vap', 'toluene')   dimensionless   0.038450   0.0088437 
Total Mole Fraction ('Vap', 'hydrogen')  dimensionless    0.16148     0.18592 
Total Mole Fraction ('Vap', 'methane')   dimensionless    0.62599     0.72074 
Temperature                                     kelvin     771.86      325.00 
Pressure                                        pascal 3.5000e+05  3.5000e+05 

8 Optimization#

We saw from the results above that the total operating cost for the base case was $419,122 per year. We are producing 0.142 mol/s of benzene at a purity of 82%. However, we are losing around 42% of benzene in F101 vapor outlet stream.

Let us try to minimize this cost such that:

  • we are producing at least 0.15 mol/s of benzene in F102 vapor outlet i.e. our product stream

  • purity of benzene i.e. the mole fraction of benzene in F102 vapor outlet is at least 80%

  • restricting the benzene loss in F101 vapor outlet to less than 20%

For this problem, our decision variables are as follows:

  • H101 outlet temperature

  • R101 cooling duty provided

  • F101 outlet temperature

  • F102 outlet temperature

  • F102 deltaP in the flash tank

Let us declare our objective function for this problem.

m.fs.objective = Objective(expr=m.fs.operating_cost)

Now, we need to unfix the decision variables as we had solved a square problem (degrees of freedom = 0) until now.

m.fs.H101.outlet.temperature.unfix()
m.fs.R101.heat_duty.unfix()
m.fs.F101.vap_outlet.temperature.unfix()
m.fs.F102.vap_outlet.temperature.unfix()
# Todo: Unfix deltaP for F102
m.fs.F102.deltaP.unfix()

Next, we need to set bounds on these decision variables to values shown below:

  • H101 outlet temperature [500, 600] K

  • R101 outlet temperature [600, 800] K

  • F101 outlet temperature [298, 450] K

  • F102 outlet temperature [298, 450] K

  • F102 outlet pressure [105000, 110000] Pa

Let us first set the variable bound for the H101 outlet temperature as shown below:

m.fs.H101.outlet.temperature[0].setlb(500)
m.fs.H101.outlet.temperature[0].setub(600)
# Todo: Set the bounds for reactor outlet temperature
m.fs.R101.outlet.temperature[0].setlb(600)
m.fs.R101.outlet.temperature[0].setub(800)

Let us fix the bounds for the rest of the decision variables.

m.fs.F101.vap_outlet.temperature[0].setlb(298.0)
m.fs.F101.vap_outlet.temperature[0].setub(450.0)
m.fs.F102.vap_outlet.temperature[0].setlb(298.0)
m.fs.F102.vap_outlet.temperature[0].setub(450.0)
m.fs.F102.vap_outlet.pressure[0].setlb(105000)
m.fs.F102.vap_outlet.pressure[0].setub(110000)

Now, the only things left to define are our constraints on overhead loss in F101, product flow rate and purity in F102. Let us first look at defining a constraint for the overhead loss in F101 where we are restricting the benzene leaving the vapor stream to less than 20 % of the benzene available in the reactor outlet.

m.fs.overhead_loss = Constraint(
    expr=m.fs.F101.control_volume.properties_out[0].flow_mol_phase_comp[
        "Vap", "benzene"
    ]
    <= 0.20
    * m.fs.R101.control_volume.properties_out[0].flow_mol_phase_comp["Vap", "benzene"]
)
# Todo: Add minimum product flow constraint
m.fs.product_flow = Constraint(
    expr=m.fs.F102.control_volume.properties_out[0].flow_mol_phase_comp[
        "Vap", "benzene"
    ]
    >= 0.15
)

Let us add the final constraint on product purity or the mole fraction of benzene in the product stream such that it is at least greater than 80%.

m.fs.product_purity = Constraint(expr=m.fs.purity >= 0.80)

We have now defined the optimization problem and we are now ready to solve this problem.

results = solver.solve(m, tee=False)

8.1 Optimization Results#

Display the results and product specifications

print("operating cost = $", value(m.fs.operating_cost))

print()
print("Product flow rate and purity in F102")

m.fs.F102.report()

print()
print("benzene purity = ", value(m.fs.purity))

print()
print("Overhead loss in F101")
m.fs.F101.report()
operating cost = $ 312674.2367537996

Product flow rate and purity in F102

====================================================================================
Unit : fs.F102                                                             Time: 0.0
------------------------------------------------------------------------------------
    Unit Performance

    Variables: 

    Key             : Value       : Units  : Fixed : Bounds
          Heat Duty :      8370.2 :   watt : False : (None, None)
    Pressure Change : -2.4500e+05 : pascal : False : (None, None)

------------------------------------------------------------------------------------
    Stream Table
                                                 Units         Inlet    Vapor Outlet  Liquid Outlet
    Total Molar Flowrate Liq                  mole / second    0.28812            -             -  
    Total Molar Flowrate Vap                  mole / second 1.0000e-08            -             -  
    Total Mole Fraction ('Liq', 'benzene')    dimensionless    0.75463            -             -  
    Total Mole Fraction ('Liq', 'toluene')    dimensionless    0.24537            -             -  
    Total Mole Fraction ('Vap', 'benzene')    dimensionless   0.032748            -             -  
    Total Mole Fraction ('Vap', 'toluene')    dimensionless  0.0032478            -             -  
    Total Mole Fraction ('Vap', 'hydrogen')   dimensionless    0.21614            -             -  
    Total Mole Fraction ('Vap', 'methane')    dimensionless    0.74786            -             -  
    Temperature                                      kelvin     301.88            -             -  
    Pressure                                         pascal 3.5000e+05            -             -  
    flow_mol_phase Liq                        mole / second          -   1.0000e-08       0.10493  
    flow_mol_phase Vap                        mole / second          -      0.18319    1.0000e-08  
    mole_frac_phase_comp ('Liq', 'benzene')   dimensionless          -      0.64256       0.64256  
    mole_frac_phase_comp ('Liq', 'toluene')   dimensionless          -      0.35744       0.35744  
    mole_frac_phase_comp ('Vap', 'benzene')   dimensionless          -      0.81883       0.81883  
    mole_frac_phase_comp ('Vap', 'toluene')   dimensionless          -      0.18117       0.18117  
    mole_frac_phase_comp ('Vap', 'hydrogen')  dimensionless          -   1.1799e-08    1.1799e-08  
    mole_frac_phase_comp ('Vap', 'methane')   dimensionless          -   4.0825e-08    4.0825e-08  
    temperature                                      kelvin          -       362.93        362.93  
    pressure                                         pascal          -   1.0500e+05    1.0500e+05  
====================================================================================

benzene purity =  0.8188295888412465

Overhead loss in F101

====================================================================================
Unit : fs.F101                                                             Time: 0.0
------------------------------------------------------------------------------------
    Unit Performance

    Variables: 

    Key             : Value   : Units  : Fixed : Bounds
          Heat Duty : -68347. :   watt : False : (None, None)
    Pressure Change :  0.0000 : pascal :  True : (None, None)

------------------------------------------------------------------------------------
    Stream Table
                                                 Units         Inlet    Vapor Outlet  Liquid Outlet
    Total Molar Flowrate Liq                  mole / second 1.5819e-12            -             -  
    Total Molar Flowrate Vap                  mole / second     1.9480            -             -  
    Total Mole Fraction ('Liq', 'benzene')    dimensionless    0.57381            -             -  
    Total Mole Fraction ('Liq', 'toluene')    dimensionless    0.42619            -             -  
    Total Mole Fraction ('Vap', 'benzene')    dimensionless    0.13952            -             -  
    Total Mole Fraction ('Vap', 'toluene')    dimensionless   0.039059            -             -  
    Total Mole Fraction ('Vap', 'hydrogen')   dimensionless    0.18417            -             -  
    Total Mole Fraction ('Vap', 'methane')    dimensionless    0.63725            -             -  
    Temperature                                      kelvin     775.95            -             -  
    Pressure                                         pascal 3.5000e+05            -             -  
    flow_mol_phase Liq                        mole / second          -   1.0000e-08       0.28812  
    flow_mol_phase Vap                        mole / second          -       1.6598    1.0000e-08  
    mole_frac_phase_comp ('Liq', 'benzene')   dimensionless          -      0.75463       0.75463  
    mole_frac_phase_comp ('Liq', 'toluene')   dimensionless          -      0.24537       0.24537  
    mole_frac_phase_comp ('Vap', 'benzene')   dimensionless          -     0.032748      0.032748  
    mole_frac_phase_comp ('Vap', 'toluene')   dimensionless          -    0.0032478     0.0032478  
    mole_frac_phase_comp ('Vap', 'hydrogen')  dimensionless          -      0.21614       0.21614  
    mole_frac_phase_comp ('Vap', 'methane')   dimensionless          -      0.74786       0.74786  
    temperature                                      kelvin          -       301.88        301.88  
    pressure                                         pascal          -   3.5000e+05    3.5000e+05  
====================================================================================

Display optimal values for the decision variables

print(
    f"""Optimal Values:

H101 outlet temperature = {value(m.fs.H101.outlet.temperature[0]):.3f} K

R101 outlet temperature = {value(m.fs.R101.outlet.temperature[0]):.3f} K

F101 outlet temperature = {value(m.fs.F101.vap_outlet.temperature[0]):.3f} K

F102 outlet temperature = {value(m.fs.F102.vap_outlet.temperature[0]):.3f} K
F102 outlet pressure = {value(m.fs.F102.vap_outlet.pressure[0]):.3f} Pa
"""
)
Optimal Values:

H101 outlet temperature = 500.000 K

R101 outlet temperature = 775.947 K

F101 outlet temperature = 301.881 K

F102 outlet temperature = 362.935 K
F102 outlet pressure = 105000.000 Pa