Skip to content

Commit 1ce6434

Browse files
committed
[Examples] Revise preconditioned integration example
- Modify formatting for Sphinx Gallery - Use a larger mechanism file from the example_data repo to better demonstrate the benefits of preconditioning
1 parent 027f224 commit 1ce6434

File tree

1 file changed

+53
-47
lines changed

1 file changed

+53
-47
lines changed

samples/python/reactors/preconditioned_integration.py

Lines changed: 53 additions & 47 deletions
Original file line numberDiff line numberDiff line change
@@ -3,7 +3,7 @@
33
========================================================================
44
55
Ideal gas, constant-pressure, adiabatic kinetics simulation that compares preconditioned
6-
and non-preconditioned integration of nDodecane.
6+
and non-preconditioned integration of n-hexane.
77
88
Requires: cantera >= 3.0.0, matplotlib >= 2.0
99
@@ -12,26 +12,26 @@
1212
import cantera as ct
1313
import numpy as np
1414
import matplotlib.pyplot as plt
15+
plt.rcParams['figure.constrained_layout.use'] = True
1516
from timeit import default_timer
1617

17-
def integrate_reactor(n_reactors=1, preconditioner=True):
18-
"""
19-
Compare the integrations of a preconditioned reactor network and a
20-
non-preconditioned reactor network. Increase the number of reactors in the
21-
network with the keyword argument `n_reactors` to see greater differences in
22-
performance.
23-
"""
24-
# Use a reduced n-dodecane mechanism with PAH formation pathways
25-
gas = ct.Solution('nDodecane_Reitz.yaml', 'nDodecane_IG')
26-
# Create multiple reactors and set initial contents
18+
# %%
19+
# Simulation setup
20+
# ----------------
21+
#
22+
# Create a reactor network for simulating the constant pressure ignition of a
23+
# stoichiometric n-hexane/air mixture, with or without the use of the preconditioned
24+
# solver.
25+
def integrate_reactor(preconditioner=True):
26+
# Use a detailed n-hexane mechanism with 1268 species
27+
gas = ct.Solution('example_data/n-hexane-NUIG-2015.yaml')
2728
gas.TP = 1000, ct.one_atm
28-
gas.set_equivalence_ratio(1, 'c12h26', 'n2:3.76, o2:1.0')
29-
reactors = [ct.IdealGasConstPressureMoleReactor(gas) for n in range(n_reactors)]
29+
gas.set_equivalence_ratio(1, 'NC6H14', 'N2:3.76, O2:1.0')
30+
reactor = ct.IdealGasConstPressureMoleReactor(gas)
3031
# set volume for reactors
31-
for r in reactors:
32-
r.volume = 0.1
32+
reactor.volume = 0.1
3333
# Create reactor network
34-
sim = ct.ReactorNet(reactors)
34+
sim = ct.ReactorNet([reactor])
3535
# Add preconditioner
3636
if preconditioner:
3737
sim.derivative_settings = {"skip-third-bodies":True, "skip-falloff":True}
@@ -40,10 +40,10 @@ def integrate_reactor(n_reactors=1, preconditioner=True):
4040
# Advance to steady state
4141
integ_time = default_timer()
4242
# solution array for state data
43-
states = ct.SolutionArray(reactors[0].thermo, extra=['time'])
43+
states = ct.SolutionArray(reactor.thermo, extra=['time'])
4444
# advance to steady state manually
4545
while (sim.time < 0.1):
46-
states.append(reactors[0].thermo.state, time=sim.time)
46+
states.append(reactor.thermo.state, time=sim.time)
4747
sim.step()
4848
integ_time = default_timer() - integ_time
4949
# Return time to integrate
@@ -53,35 +53,41 @@ def integrate_reactor(n_reactors=1, preconditioner=True):
5353
print(f"Non-preconditioned Integration Time: {integ_time:f}")
5454
# Get and output solver stats
5555
for key, value in sim.solver_stats.items():
56-
print(key, value)
56+
print(f"{key:>24s}: {value}")
5757
print("\n")
5858
# return some variables for plotting
59-
return states.time, states.T, states('CO2').Y, states('C12H26').Y
59+
return states.time, states.T, states('CO2').Y, states('NC6H14').Y
6060

61-
if __name__ == "__main__":
62-
# integrate both to steady state
63-
timep, Tp, CO2p, C12H26p = integrate_reactor(preconditioner=True)
64-
timenp, Tnp, CO2np, C12H26np = integrate_reactor(preconditioner=False)
65-
# plot selected state variables
66-
fig, axs = plt.subplots(1, 3)
67-
fig.tight_layout()
68-
ax1, ax2, ax3 = axs
69-
# temperature plot
70-
ax1.set_xlabel("Time")
71-
ax1.set_ylabel("Temperature")
72-
ax1.plot(timenp, Tnp, linewidth=2)
73-
ax1.plot(timep, Tp, linewidth=2, linestyle=":")
74-
ax1.legend(["Normal", "Preconditioned"])
75-
# CO2 plot
76-
ax2.set_xlabel("Time")
77-
ax2.set_ylabel("CO2")
78-
ax2.plot(timenp, CO2np, linewidth=2)
79-
ax2.plot(timep, CO2p, linewidth=2, linestyle=":")
80-
ax2.legend(["Normal", "Preconditioned"])
81-
# C12H26 plot
82-
ax3.set_xlabel("Time")
83-
ax3.set_ylabel("C12H26")
84-
ax3.plot(timenp, C12H26np, linewidth=2)
85-
ax3.plot(timep, C12H26p, linewidth=2, linestyle=":")
86-
ax3.legend(["Normal", "Preconditioned"])
87-
plt.show()
61+
# %%
62+
# Integrate with sparse, preconditioned solver
63+
# --------------------------------------------
64+
timep, Tp, CO2p, NC6H14p = integrate_reactor(preconditioner=True)
65+
66+
# %%
67+
# Integrate with direct linear solver
68+
# -----------------------------------
69+
timenp, Tnp, CO2np, NC6H14np = integrate_reactor(preconditioner=False)
70+
71+
# %%
72+
# Plot selected state variables
73+
# -----------------------------
74+
fig, (ax1, ax2, ax3) = plt.subplots(3, 1, figsize=(5, 8))
75+
# temperature plot
76+
ax1.set_xlabel("Time")
77+
ax1.set_ylabel("Temperature")
78+
ax1.plot(timenp, Tnp, linewidth=2)
79+
ax1.plot(timep, Tp, linewidth=2, linestyle=":")
80+
ax1.legend(["Normal", "Preconditioned"])
81+
# CO2 plot
82+
ax2.set_xlabel("Time")
83+
ax2.set_ylabel("CO2")
84+
ax2.plot(timenp, CO2np, linewidth=2)
85+
ax2.plot(timep, CO2p, linewidth=2, linestyle=":")
86+
ax2.legend(["Normal", "Preconditioned"])
87+
# C12H26 plot
88+
ax3.set_xlabel("Time")
89+
ax3.set_ylabel("NC6H14")
90+
ax3.plot(timenp, NC6H14np, linewidth=2)
91+
ax3.plot(timep, NC6H14p, linewidth=2, linestyle=":")
92+
ax3.legend(["Normal", "Preconditioned"])
93+
plt.show()

0 commit comments

Comments
 (0)