Figure 7.32:

Pellet CO profiles at several reactor positions.

Figure 7.32

Code for Figure 7.32

Text of the GNU GPL.

main.py


  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
# Converted from multicolloc_log.m - multi-component packed bed with pellet collocation
# Uses SUNDIALS IDA (via scikits.odes) to solve the DAE directly, matching ode15i
import numpy as np
from misc import colloc, octave_save
from scikits.odes import dae

Rg   = 8.314; Mair = 28.96
Pin  = 2.0*1.013e5; Pext = 1.013e5
Tin  = 550.
Rt   = 5.0; At = np.pi*Rt**2
Ta   = 325.; U = 5.5e-3
delH1 = -67.63e3; delH2 = -460.4e3
Cp   = 0.25; vis = 0.028e-2
bt   = 1.; u = 75./bt
Ac   = bt*At; Qin = u*Ac
cfin = Pin/(Rg*Tin)*1e-6
Rp   = 0.175; a = Rp/3.
rhob = 0.51; rhop = 0.68
epsb = 1 - rhob/rhop
xc   = 0.996

cafin = cfin*0.02; cbfin = cfin*0.03; ccfin = cfin*5e-4
Nafin = Qin*cafin; Nbfin = Qin*cbfin; Ncfin = Qin*ccfin
Nfin  = np.array([Nafin, Nbfin, Ncfin])
Ntfin = Qin*cfin
massf = Ntfin*Mair
Cptot = massf*Cp

k100 = 6.802e16*2.6e6*80./100*0.05/100
k200 = 1.416e18*2.6e6*80./100*0.05/100
E1=13108.; E2=15109.
Ka0=8.099e6; Kc0=2.579e8
Ea=-409.; Ec=191.
Da=0.0487; Db=0.0469; Dc=0.0487
kma=0.4*9.76; kmb=0.4*10.18; kmc=0.4*9.76

npts = 40
R_col, A_col, B_col, Q_col = colloc(npts-2, 'left', 'right')
R_col = R_col*Rp; A_col = A_col/Rp; B_col = B_col/Rp**2

p = dict(npts=npts, A=A_col, B=B_col, R=R_col,
         k10=k100, k20=k200, E1=E1, E2=E2,
         Ea=Ea, Ec=Ec, Ka0=Ka0, Kc0=Kc0,
         Da=Da, Db=Db, Dc=Dc, kma=kma, kmb=kmb, kmc=kmc,
         T=Tin, Nafin=Nafin, Ncfin=Ncfin, Ntfin=Ntfin,
         Pin=Pin, Tin=Tin, epsb=epsb, a=a, Qin=Qin,
         Cptot=Cptot, U=U, Ta=Ta, Rt=Rt,
         delH1=delH1, delH2=delH2, Rp=Rp, vis=vis, massf=massf, Ac=Ac,
         xc=xc, Pext=Pext, caf=cafin, cbf=cbfin, ccf=ccfin)

def pelletfull(x, p):
    n  = p['npts']
    za = x[:n]; zb = x[n:2*n]; zc = x[2*n:3*n]
    ca = np.exp(za); cb = np.exp(zb); cc = np.exp(zc)
    k1 = p['k10']*np.exp(-p['E1']/p['T'])
    k2 = p['k20']*np.exp(-p['E2']/p['T'])
    Ka = p['Ka0']*np.exp(-p['Ea']/p['T'])
    Kc = p['Kc0']*np.exp(-p['Ec']/p['T'])
    den = (1+Ka*ca+Kc*cc)**2
    r1 = k1*ca*cb/den; r2 = k2*cc*cb/den
    Ra = -r1; Rb = -0.5*r1 - 4.5*r2; Rc_ = -r2
    resid = np.zeros(3*n)
    resid[0:n] = p['B']@za + (p['A']@za)*(p['A']@za + 2./p['R']) + Ra/(p['Da']*ca)
    resid[0]   = p['A'][0,:]@za
    resid[n-1] = p['Da']*(p['A'][-1,:]@za) - p['kma']*(p['caf']/ca[-1] - 1)
    resid[n:2*n] = p['B']@zb + (p['A']@zb)*(p['A']@zb + 2./p['R']) + Rb/(p['Db']*cb)
    resid[n]     = p['A'][0,:]@zb
    resid[2*n-1] = p['Db']*(p['A'][-1,:]@zb) - p['kmb']*(p['cbf']/cb[-1] - 1)
    resid[2*n:3*n] = p['B']@zc + (p['A']@zc)*(p['A']@zc + 2./p['R']) + Rc_/(p['Dc']*cc)
    resid[2*n]     = p['A'][0,:]@zc
    resid[3*n-1]   = p['Dc']*(p['A'][-1,:]@zc) - p['kmc']*(p['ccf']/cc[-1] - 1)
    rate = np.array([p['A'][-1,:]@ca, p['A'][-1,:]@cb, p['A'][-1,:]@cc])
    return resid, rate

def pellet_resid(x, p):
    return pelletfull(x, p)[0]

# Homotopy to find inlet pellet profile
from scipy.optimize import fsolve
p['caf'] = cafin; p['cbf'] = cbfin; p['ccf'] = ccfin
kstep=0.1; ksuc=0.; ktry=0.; kmin=1e-6
z0 = np.log(np.concatenate([cafin*np.ones(npts), cbfin*np.ones(npts), ccfin*np.ones(npts)]))
while ksuc < 1.0:
    ktry = min(1.0, ktry+kstep)
    p['k10'] = ktry*k100; p['k20'] = ktry*k200
    z, _, info, _ = fsolve(pellet_resid, z0, args=(p,), full_output=True,
                           maxfev=50000, factor=0.1)[:4]
    if info != 1:
        ktry = ksuc; kstep /= 2
        if kstep < kmin: raise RuntimeError('homotopy failed')
    else:
        ksuc = ktry; z0 = z
p['k10'] = k100; p['k20'] = k200

# DAE residual — mirrors the Octave bed() function used with ode15i
# State:  y = [Naf, Nbf, Ncf, T, P, z_pellet...]   (5 + 3*npts states)
# ydot:   time-derivatives of all states
# res[0:5]  = ydot[0:5] - rhs_bed   (differential equations)
# res[5:]   = pelletres(z, pp)       (algebraic constraints)
nbed = 5
npel = 3 * npts

def dae_residual(V, y, ydot, res):
    Naf, Nbf, Ncf, T, P = y[:nbed]
    z_pel = y[nbed:]
    Ntf = p['Ntfin'] - 0.5*(p['Nafin']-Naf) + 0.5*(p['Ncfin']-Ncf)
    Q   = p['Qin']*p['Pin']/P * T/p['Tin'] * Ntf/p['Ntfin']
    pp = dict(p)
    pp['caf'] = Naf/Q; pp['cbf'] = Nbf/Q; pp['ccf'] = Ncf/Q
    pp['T'] = T
    pel_res, rate = pelletfull(z_pel, pp)
    dcadr, dcbdr, dccdr = rate
    r1p = Da/a*dcadr; r2p = Dc/a*dccdr
    Rj  = -(1-epsb)/a * np.array([Da*dcadr, Db*dcbdr, Dc*dccdr])
    dT  = (-(delH1*r1p+delH2*r2p) + 2./Rt*U*(Ta-T)) / Cptot
    dP  = (-(1-epsb)*Q / (2*Rp*epsb**3*Ac**2) *
           (150*(1-epsb)*vis/(2*Rp) + 7./4*massf/Ac))
    f_bed = np.concatenate([Rj, [dT, dP]])
    # differential residuals
    res[:nbed] = ydot[:nbed] - f_bed
    # algebraic residuals (pellet BVP)
    res[nbed:] = pel_res

nvs=201; vfinal=2000.
vsteps = np.linspace(0., vfinal, nvs)
y0    = np.concatenate([np.array([Nafin, Nbfin, Ncfin, Tin, Pin]), z0])
ydot0 = np.zeros(nbed + npel)

# algebraic components: indices nbed .. nbed+npel-1
# IDA needs to know which are differential (1) vs algebraic (0)
algvar = np.ones(nbed + npel)
algvar[nbed:] = 0.

solver = dae('ida', dae_residual,
             algebraic_vars_idx=np.where(algvar == 0)[0],
             compute_initcond='yp0',   # let IDA compute consistent ydot0
             old_api=False,
             rtol=1e-7, atol=1e-7)

solution = solver.solve(vsteps, y0, ydot0)

# Check for termination event (xc conversion or pressure drop)
sol_t = solution.values.t
sol_y = solution.values.y   # shape (nout, nstates)

# Trim at stopping condition
mask = np.ones(len(sol_t), dtype=bool)
for i, (t, y) in enumerate(zip(sol_t, sol_y)):
    conv = 1 - y[2]/p['Ncfin']
    if conv >= p['xc'] or y[4] <= p['Pext']:
        mask[i+1:] = False
        break
sol_t = sol_t[mask]
sol_y = sol_y[mask]

nout = len(sol_t); vout = sol_t
y_out = sol_y[:, :nbed]
P_out=y_out[:,4]; T_out=y_out[:,3]
ctot = P_out/(Rg*T_out)*1e-6
Nt_out = Ntfin - 0.5*(Nafin-y_out[:,0]) + 0.5*(Ncfin-y_out[:,2])
Q_out  = Nt_out/ctot
caf_out=y_out[:,0]/Q_out; cbf_out=y_out[:,1]/Q_out; ccf_out=y_out[:,2]/Q_out
Patm_out = P_out/1.013e5
table1 = np.column_stack([vout, caf_out, cbf_out, ccf_out, T_out, Patm_out])

# Pellet profiles at selected rows — read directly from state vector
rowsp = np.array([0, 4, 9, 49, 89, 99, nout-1])
rowsp = rowsp[rowsp < nout]
ca_cols=[]; cb_cols=[]; cc_cols=[]
for ri in rowsp:
    z_ri = sol_y[ri, nbed:]
    ca_cols.append(np.exp(z_ri[:npts]).reshape(-1,1))
    cb_cols.append(np.exp(z_ri[npts:2*npts]).reshape(-1,1))
    cc_cols.append(np.exp(z_ri[2*npts:3*npts]).reshape(-1,1))

table2 = np.column_stack([R_col] + ca_cols + cb_cols + cc_cols)
octave_save('multicolloc_log.dat', ('table1', table1), ('table2', table2))