Figure 6.24:

Real and imaginary parts of the eigenvalues of the Jacobian matrix near points C--F.

Figure 6.24

Code for Figure 6.24

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
# Converted from st_st_osc_eigo.m - Oscillatory system eigenvalues (upper/limit-cycle branch)
# Saves rows 441:rowend (the upper/limit-cycle segment)
import numpy as np
from scipy.optimize import fsolve
from misc import save_ascii

p = dict(
    k_m      = 0.004,
    T_m      = 298.,
    E        = 15000.,
    c_Af     = 2.,
    C_p      = 4.,
    rho      = 1000.,
    T_f      = 298.,
    T_a      = 298.,
    DeltaH_R = -2.2e5,
    U        = 340.,
    theta    = 35.,
    T_set    = 321.53,
    T_fs     = 298.,
    Kc       = 0.,
)
p['C_ps'] = p['C_p'] * p['rho']

def st_st_cA(x, c_A, p):
    theta = x[0]; T = x[1]
    k = p['k_m'] * np.exp(-p['E']*(1./T - 1./p['T_m']))
    return [p['c_Af'] - (1+k*theta)*c_A,
            p['U']*theta*(p['T_a']-T) + p['C_ps']*(p['T_f']-T) - k*theta*c_A*p['DeltaH_R']]

def st_st_theta(x, theta, p):
    c_A = x[0]; T = x[1]
    k = p['k_m'] * np.exp(-p['E']*(1./T - 1./p['T_m']))
    return [p['c_Af'] - (1+k*theta)*c_A,
            p['U']*theta*(p['T_a']-T) + p['C_ps']*(p['T_f']-T) - k*theta*c_A*p['DeltaH_R']]

def compute_lamrow(k, c_A, T, theta, p):
    Jac = np.array([[-1./theta - k, -k*c_A*p['E']/(T*T)],
                    [-k*p['DeltaH_R']/p['C_ps'],
                     -p['U']/p['C_ps'] - 1./theta - k*c_A*p['DeltaH_R']/p['C_ps']*p['E']/(T*T)]])
    lam = np.linalg.eigvals(Jac)
    a = lam[0].real; b = lam[0].imag
    c = lam[1].real; d = lam[1].imag
    return [a,b,c,d] if a >= c else [c,d,a,b]

# Lower branch (logspace in conversion)
npts1 = 200
x0 = [1., p['T_f']]
xvect = np.logspace(np.log10(1e-4), np.log10(0.75), npts1)
c_Avect = (1 - xvect) * p['c_Af']
table = []

for c_A in c_Avect:
    X, _, ier, _ = fsolve(lambda x: st_st_cA(x, c_A, p), x0, full_output=True)
    theta = X[0]; T = X[1]
    conv  = (p['c_Af'] - c_A) / p['c_Af']
    k = p['k_m'] * np.exp(-p['E']*(1./T - 1./p['T_m']))
    lamrow = compute_lamrow(k, c_A, T, theta, p)
    table.append([theta, T, conv] + lamrow + [float(ier)])
    x0 = X

c_A_c = c_A; T_c = T; theta_c = theta; conv_c = conv

# Near-extinction patch
npts2 = 50
xvect2 = np.linspace(conv_c, 0.95, npts2)
c_Avect2 = (1 - xvect2) * p['c_Af']

for c_A in c_Avect2:
    X, _, ier, _ = fsolve(lambda x: st_st_cA(x, c_A, p), x0, full_output=True)
    theta = X[0]; T = X[1]
    conv  = (p['c_Af'] - c_A) / p['c_Af']
    k = p['k_m'] * np.exp(-p['E']*(1./T - 1./p['T_m']))
    lamrow = compute_lamrow(k, c_A, T, theta, p)
    table.append([theta, T, conv] + lamrow + [float(ier)])
    x0 = X

c_A_c2 = c_A; T_c2 = T; theta_c2 = theta

# Upper branch segment 1 (logspace theta from corner to 20.7)
npts3 = 200
x0 = [c_A_c2, T_c2]
theta_vect3 = np.logspace(np.log10(theta_c2), np.log10(20.7), npts3)

for theta in theta_vect3:
    X, _, ier, _ = fsolve(lambda x: st_st_theta(x, theta, p), x0, full_output=True)
    c_A = X[0]; T = X[1]
    conv = (p['c_Af'] - c_A) / p['c_Af']
    k = p['k_m'] * np.exp(-p['E']*(1./T - 1./p['T_m']))
    lamrow = compute_lamrow(k, c_A, T, theta, p)
    table.append([theta, T, conv] + lamrow + [float(ier)])
    x0 = X

theta_c3 = theta; c_A_c3 = X[0]; T_c3 = X[1]

# Fine patch near theta=[20.7, 20.8]
npts4 = 100
x0 = [c_A_c3, T_c3]
theta_vect4 = np.linspace(theta_c3, 20.8, npts4)

for theta in theta_vect4:
    X, _, ier, _ = fsolve(lambda x: st_st_theta(x, theta, p), x0, full_output=True)
    c_A = X[0]; T = X[1]
    conv = (p['c_Af'] - c_A) / p['c_Af']
    k = p['k_m'] * np.exp(-p['E']*(1./T - 1./p['T_m']))
    lamrow = compute_lamrow(k, c_A, T, theta, p)
    table.append([theta, T, conv] + lamrow + [float(ier)])
    x0 = X

theta_c4 = theta; c_A_c4 = X[0]; T_c4 = X[1]

# Finish to large theta
npts5 = 200
x0 = [c_A_c4, T_c4]
theta_vect5 = np.logspace(np.log10(theta_c4), np.log10(500.), npts5)

for theta in theta_vect5:
    X, _, ier, _ = fsolve(lambda x: st_st_theta(x, theta, p), x0, full_output=True)
    c_A = X[0]; T = X[1]
    conv = (p['c_Af'] - c_A) / p['c_Af']
    k = p['k_m'] * np.exp(-p['E']*(1./T - 1./p['T_m']))
    lamrow = compute_lamrow(k, c_A, T, theta, p)
    table.append([theta, T, conv] + lamrow + [float(ier)])
    x0 = X

table = np.array(table)
rowend = npts1 + npts2 + npts3 + npts4 + npts5

# Octave saves rows 441:rowend (1-indexed), i.e. Python 440:rowend
tmp = table[440:rowend, :]
save_ascii('st_st_osc_eigo.dat', tmp)