Ignis
Loading...
Searching...
No Matches
c2h4red.cc
Go to the documentation of this file.
1
13#include <cmath>
14#include <algorithm>
15
16using namespace std;
17
18/* Table of constant values */
19
20static double c_b5 = 10.;
21
22
23/* 15-step Reduced Mechanism for C2H4/Air */
24
25/* August 25, 2005 */
26
27/* Developed by Tianfeng Lu */
28/* MAE Dept., Princeton Univ. */
29/* D103 E-Quad, Olden St */
30/* Princeton, NJ 08544 */
31/* Email: tlu@princeton.edu */
32
33/* Applicable parameter range: */
34/* Equivalence ratio: 0.5-1.7 (for premixed mixtures) */
35/* Atmospheric or slightly higer pressure */
36/* Ingition, extinction, and flames */
37/* C */
38/* ----------------------------------------------------------------------C */
39
41// int *ickwrk, double *rckwrk, double *rrsp)
42
43/* Subroutine */
44
45// int getrates_(double *p, double *t, double *y, double *rrsp)
46
47void getProblemSpecificRR(double rho, double temp, double pres,
48 double *yi, double *rrsp) {
49
50 // rho not used but fits the generic interface: units kg/m3
51 // temp units are K
52 // rrsp is output: units are kmol/m^3*s
53 // pres input units are Pa
54
55 pres *= 10.0; // convert Pa to dynes/cm^2
56
57
58 /* Initialized data */
59
60 static double small = 1e-50;
61 static double rf[167] = { 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
62 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
63 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
64 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
65 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
66 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
67 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
68 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0. };
69 static double rb[167] = { 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
70 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
71 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
72 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
73 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
74 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
75 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,
76 0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0.,0. };
77 static double xq[10] = { 0.,0.,0.,0.,0.,0.,0.,0.,0.,0. };
78
79 /* System generated locals */
80 double d__1, d__2, d__3;
81
82 /* Builtin functions */
83 //dol double log(double);
84
85 /* Local variables */
86 static double c__[19];
87 static int i__, k;
88 static double a1, a2, a3, a4, a5, a6, a9, b10, a12, a21, a23, a24,
89 a25, a32, a34, a35, a29, a42, a43, a45, a46, a49, a52, a53, a54,
90 a64, a92, a94, b12, b13, b14, b20, b22, b23, b24, b30, b32, b33,
91 b34, c12, c13, c10, c22, c23, c20, den[10], abv[10], rop[167],
92 sum;
93 extern /* Subroutine */ int ratt_(double *, double *, double *
94 , double *, double *), ratx_(double *, double *,
95 double *, double *, double *, double *);
96 static double alogt, rklow[18];
97
98
99 /* Parameter adjustments */
100 --rrsp;
101// --rckwrk;
102// --ickwrk;
103 --yi;
104
105 /* Function Body */
106
107/* Convert mass fraction to mole concentration */
108
109 c__[0] = yi[1] * .496046521;
110 c__[1] = yi[2] * .992093043;
111 c__[2] = yi[3] * .0625023433;
112 c__[3] = yi[4] * .0312511716;
113 c__[4] = yi[5] * .0587980383;
114 c__[5] = yi[6] * .0555082499;
115 c__[6] = yi[7] * .0302968146;
116 c__[7] = yi[8] * .0293990192;
117 c__[8] = yi[9] * .0665112065;
118 c__[9] = yi[10] * .0623323639;
119 c__[10] = yi[11] * .0357008335;
120 c__[11] = yi[12] * .0227221341;
121 c__[12] = yi[13] * .0333039255;
122 c__[13] = yi[14] * .0384050525;
123 c__[14] = yi[15] * .0356453112;
124 c__[15] = yi[16] * .0332556033;
125 c__[16] = yi[17] * .0237882046;
126 c__[17] = yi[18] * .0237635408;
127 c__[18] = yi[19] * .0356972032;
128
129 sum = 0.f;
130 for (k = 1; k <= 19; ++k) {
131 sum += c__[k - 1];
132 }
133 sum = pres / (sum * temp * 83145100.);
134
135 for (k = 1; k <= 19; ++k) {
136/* Computing MAX */
137 d__1 = c__[k - 1];
138 c__[k - 1] = max(d__1,small);
139 c__[k - 1] *= sum;
140 }
141
142/* Compute T- and P- dependent elementary reaction rates */
143
144 alogt = log(temp);
145 ratt_(&temp, &alogt, rf, rb, rklow);
146 ratx_(&temp, &alogt, c__, rf, rb, rklow);
147
148/* Compute QSS species concetrations */
149
150/* C */
151 abv[0] = rb[72] + rb[97] + rb[98];
152 den[0] = rb[42] + rf[72] + rf[97] + rf[98];
153/* CH */
154 abv[1] = rb[5] + rb[101] + rb[102] + rb[104] + rb[107] + rb[108] + rb[149]
155 ;
156 den[1] = rf[5] + rf[42] + rb[44] + rf[73] + rb[75] + rf[99] + rf[100] +
157 rf[101] + rf[104] + rf[105] + rf[106] + rf[107] + rf[149];
158/* CH2 */
159 abv[2] = rf[17] + rf[24] + rb[43] + rb[74] + rf[77] + rb[92] + rb[102] +
160 rb[109] + rb[110] + rb[111] + rb[111] + rb[112] + rb[113] + rb[
161 114] + rb[150] + rb[151] + rb[152] + rb[152];
162 den[2] = rf[6] + rb[17] + rb[24] + rf[43] + rf[74] + rf[75] + rb[77] + rf[
163 92] + rb[100] + rf[109] + rf[110] + rf[112] + rf[113] + rf[114] +
164 rb[116] + rb[120] + rb[123] + rb[124] + rf[150] + rf[151];
165/* CH2 */
166 abv[3] = rb[7] + rb[76] + rf[78] + rb[117] + rb[118] + rb[119] + rb[121]
167 + rb[122] + rb[125] + rb[153];
168 den[3] = rf[7] + rf[8] + rf[44] + rb[53] + rb[62] + rf[76] + rb[78] + rf[
169 116] + rf[117] + rf[118] + rf[119] + rf[120] + rf[121] + rf[122]
170 + rf[123] + rf[124] + rf[125] + rf[153];
171/* HCO */
172 abv[4] = rb[12] + rb[13] + rf[14] + rf[19] + rf[26] + rb[47] + rb[48] +
173 rf[50] + rb[81] + rf[82] + rf[96] + rb[132] + rf[133] + rb[136] +
174 rb[137] + rb[138] + rf[165];
175 den[4] = rb[6] + rb[8] + rf[12] + rf[13] + rb[14] + rb[19] + rb[26] + rf[
176 47] + rf[48] + rb[50] + rb[73] + rf[81] + rb[82] + rb[96] + rb[99]
177 + rb[106] + rf[132] + rb[133] + rf[136] + rf[137] + rf[138] + rb[
178 160];
179/* CH3O */
180 abv[5] = rb[15] + rf[49] + rb[51] + rb[52] + rb[83] + rf[94] + rf[127] +
181 rb[139];
182 den[5] = rf[15] + rb[49] + rf[51] + rf[52] + rf[53] + rf[83] + rb[94] +
183 rb[127] + rf[139];
184/* C2H3 */
185 abv[6] = rb[18] + rf[54] + rb[55] + rb[56] + rf[58] + rb[86] + rf[87] +
186 rf[134] + rb[155] + rb[166];
187 den[6] = rf[18] + rb[54] + rf[55] + rf[56] + rb[58] + rf[86] + rb[87] +
188 rb[134] + rf[140] + rf[154] + rf[155] + rf[166];
189/* C2H5 */
190 abv[7] = rb[20] + rf[21] + rf[57] + rb[59] + rb[60] + rf[61] + rf[88] +
191 rf[131] + rf[135] + rb[142] + rf[165];
192 den[7] = rf[20] + rb[21] + rb[57] + rf[59] + rf[60] + rb[61] + rb[88] +
193 rb[131] + rb[135] + rf[142];
194/* HCCO */
195 abv[8] = rf[16] + rb[22] + rf[23] + rf[63] + rf[89] + rb[108] + rb[143] +
196 rb[144] + rb[144];
197 den[8] = rb[16] + rf[22] + rb[23] + rf[62] + rb[63] + rb[89] + rb[105] +
198 rf[143];
199/* CH2CHO */
200 abv[9] = rf[146] + rf[156] + rb[158] + rb[161] + rb[162];
201 den[9] = rb[146] + rb[156] + rf[158] + rf[161] + rf[162];
202
203/* C2H3 */
204 xq[6] = abv[6] / den[6];
205/* C2H5 */
206 xq[7] = abv[7] / den[7];
207
208 a1 = abv[0] / den[0];
209 a12 = rf[42] / den[0];
210 a2 = abv[1] / den[1];
211 a21 = rb[42] / den[1];
212 a23 = (rf[75] + rb[100]) / den[1];
213 a24 = rf[44] / den[1];
214 a25 = (rb[73] + rb[99] + rb[106]) / den[1];
215 a29 = rb[105] / den[1];
216 a3 = abv[2] / den[2];
217 a32 = (rb[75] + rf[100]) / den[2];
218 a34 = (rf[116] + rf[120] + rf[123] + rf[124]) / den[2];
219 a35 = rb[6] / den[2];
220 a4 = abv[3] / den[3];
221 a42 = rb[44] / den[3];
222 a43 = (rb[116] + rb[120] + rb[123] + rb[124]) / den[3];
223 a45 = rb[8] / den[3];
224 a46 = rf[53] / den[3];
225 a49 = rf[62] / den[3];
226 a5 = (abv[4] + rf[140] * xq[6]) / den[4];
227 a52 = (rf[73] + rf[99] + rf[106]) / den[4];
228 a53 = rf[6] / den[4];
229 a54 = rf[8] / den[4];
230 a6 = abv[5] / den[5];
231 a64 = rb[53] / den[5];
232 a9 = abv[8] / den[8];
233 a92 = rf[105] / den[8];
234 a94 = rb[62] / den[8];
235
236 b10 = a2 + a21 * a1 + a5 * a25 + a9 * a29;
237 b12 = a21 * a12 + a25 * a52 + a29 * a92 - 1;
238 b13 = a23 + a25 * a53;
239 b14 = a24 + a25 * a54 + a29 * a94;
240 b20 = a3 + a35 * a5;
241 b22 = a32 + a35 * a52;
242 b23 = a35 * a53 - 1;
243 b24 = a34 + a35 * a54;
244 b30 = a4 + a45 * a5 + a46 * a6 + a49 * a9;
245 b32 = a42 + a45 * a52 + a49 * a92;
246 b33 = a43 + a45 * a53;
247 b34 = a45 * a54 + a46 * a64 + a49 * a94 - 1;
248 c12 = b34 * b12 - b14 * b32;
249 c13 = b34 * b13 - b14 * b33;
250 c10 = b34 * b10 - b14 * b30;
251 c22 = b34 * b22 - b24 * b32;
252 c23 = b34 * b23 - b24 * b33;
253 c20 = b34 * b20 - b24 * b30;
254
255/* CH */
256/* Computing MAX */
257 d__3 = (d__2 = c23 * c12 - c13 * c22, abs(d__2));
258 xq[1] = (d__1 = c13 * c20 - c23 * c10, abs(d__1)) / max(d__3,small);
259/* CH2 */
260/* Computing MAX */
261 d__2 = abs(c23);
262 xq[2] = (d__1 = c22 * xq[1] + c20, abs(d__1)) / max(d__2,small);
263/* CH2* */
264/* Computing MAX */
265 d__2 = abs(b34);
266 xq[3] = (d__1 = b32 * xq[1] + b33 * xq[2] + b30, abs(d__1)) / max(d__2,
267 small);
268/* C */
269 xq[0] = a1 + a12 * xq[1];
270/* HCO */
271 xq[4] = a5 + a52 * xq[1] + a53 * xq[2] + a54 * xq[3];
272/* CH3O */
273 xq[5] = a6 + a64 * xq[3];
274/* HCCO */
275 xq[8] = a9 + a92 * xq[1] + a94 * xq[3];
276/* CH2CHO */
277 xq[9] = (abv[9] + rb[160] * xq[4] + rf[154] * xq[6]) / den[9];
278
279/* Multiply QSS species concetratios to involved reactions */
280
281 rf[5] *= xq[1];
282 rf[6] *= xq[2];
283 rb[6] *= xq[4];
284 rf[7] *= xq[3];
285 rf[8] *= xq[3];
286 rb[8] *= xq[4];
287 rf[12] *= xq[4];
288 rf[13] *= xq[4];
289 rb[14] *= xq[4];
290 rf[15] *= xq[5];
291 rb[16] *= xq[8];
292 rb[17] *= xq[2];
293 rf[18] *= xq[6];
294 rb[19] *= xq[4];
295 rf[20] *= xq[7];
296 rb[21] *= xq[7];
297 rf[22] *= xq[8];
298 rb[23] *= xq[8];
299 rb[24] *= xq[2];
300 rb[26] *= xq[4];
301 rf[42] *= xq[1];
302 rb[42] *= xq[0];
303 rf[43] *= xq[2];
304 rf[44] *= xq[3];
305 rb[44] *= xq[1];
306 rf[47] *= xq[4];
307 rf[48] *= xq[4];
308 rb[49] *= xq[5];
309 rb[50] *= xq[4];
310 rf[51] *= xq[5];
311 rf[52] *= xq[5];
312 rf[53] *= xq[5];
313 rb[53] *= xq[3];
314 rb[54] *= xq[6];
315 rf[55] *= xq[6];
316 rf[56] *= xq[6];
317 rb[57] *= xq[7];
318 rb[58] *= xq[6];
319 rf[59] *= xq[7];
320 rf[60] *= xq[7];
321 rb[61] *= xq[7];
322 rf[62] *= xq[8];
323 rb[62] *= xq[3];
324 rb[63] *= xq[8];
325 rf[72] *= xq[0];
326 rf[73] *= xq[1];
327 rb[73] *= xq[4];
328 rf[74] *= xq[2];
329 rf[75] *= xq[2];
330 rb[75] *= xq[1];
331 rf[76] *= xq[3];
332 rb[77] *= xq[2];
333 rb[78] *= xq[3];
334 rf[81] *= xq[4];
335 rb[82] *= xq[4];
336 rf[83] *= xq[5];
337 rf[86] *= xq[6];
338 rb[87] *= xq[6];
339 rb[88] *= xq[7];
340 rb[89] *= xq[8];
341 rf[92] *= xq[2];
342 rb[94] *= xq[5];
343 rb[96] *= xq[4];
344 rf[97] *= xq[0];
345 rf[98] *= xq[0];
346 rf[99] *= xq[1];
347 rb[99] *= xq[4];
348 rf[100] *= xq[1];
349 rb[100] *= xq[2];
350 rf[101] *= xq[1];
351 rf[104] *= xq[1];
352 rf[105] *= xq[1];
353 rb[105] *= xq[8];
354 rf[106] *= xq[1];
355 rb[106] *= xq[4];
356 rf[107] *= xq[1];
357 rf[109] *= xq[2];
358 rf[110] *= xq[2];
359 rf[112] *= xq[2];
360 rf[113] *= xq[2];
361 rf[114] *= xq[2];
362 rf[116] *= xq[3];
363 rb[116] *= xq[2];
364 rf[117] *= xq[3];
365 rf[118] *= xq[3];
366 rf[119] *= xq[3];
367 rf[120] *= xq[3];
368 rb[120] *= xq[2];
369 rf[121] *= xq[3];
370 rf[122] *= xq[3];
371 rf[123] *= xq[3];
372 rb[123] *= xq[2];
373 rf[124] *= xq[3];
374 rb[124] *= xq[2];
375 rf[125] *= xq[3];
376 rb[127] *= xq[5];
377 rb[131] *= xq[7];
378 rf[132] *= xq[4];
379 rb[133] *= xq[4];
380 rb[134] *= xq[6];
381 rb[135] *= xq[7];
382 rf[136] *= xq[4];
383 rf[137] *= xq[4];
384 rf[138] *= xq[4];
385 rf[139] *= xq[5];
386 rf[140] *= xq[6];
387 rf[142] *= xq[7];
388 rf[143] *= xq[8];
389 rb[146] *= xq[9];
390 rf[149] *= xq[1];
391 rf[150] *= xq[2];
392 rf[151] *= xq[2];
393 rf[153] *= xq[3];
394 rf[154] *= xq[6];
395 rf[155] *= xq[6];
396 rb[156] *= xq[9];
397 rf[158] *= xq[9];
398 rb[160] *= xq[4];
399 rf[161] *= xq[9];
400 rf[162] *= xq[9];
401 rf[166] *= xq[6];
402
403 for (i__ = 1; i__ <= 167; ++i__) {
404 rop[i__ - 1] = rf[i__ - 1] - rb[i__ - 1];
405 }
406
407 rop[27] += rop[28];
408 rop[27] += rop[29];
409 rop[27] += rop[30];
410 rop[32] += rop[33];
411 rop[32] += rop[34];
412 rop[32] += rop[35];
413 rop[69] += rop[147];
414 rop[70] += rop[71];
415 rop[90] += rop[91];
416 rop[116] += rop[120];
417 rop[116] += rop[123];
418 rop[116] += rop[124];
419 rop[136] += rop[137];
420
421 rrsp[1] = -rop[2] + rop[7] + rop[32] + rop[38] + rop[40] + rop[42] + rop[
422 44] + rop[46] + rop[48] + rop[50] + rop[51] + rop[56] + rop[58] +
423 rop[60] + rop[61] + rop[63] - rop[65] - rop[66] - rop[100] - rop[
424 110] + rop[111] - rop[119] + rop[141] + rop[145] + rop[148] - rop[
425 149] + rop[153] + rop[161];
426 rrsp[2] = -rop[1] + rop[2] + rop[5] + rop[6] + rop[8] + rop[9] + rop[13]
427 + rop[16] + rop[18] + rop[22] - rop[27] - rop[31] - rop[32] - rop[
428 32] - rop[36] - rop[37] - rop[38] - rop[39] - rop[40] - rop[41] -
429 rop[42] - rop[43] - rop[44] - rop[45] - rop[46] - rop[47] - rop[
430 48] - rop[49] - rop[50] - rop[51] - rop[52] - rop[53] - rop[54] -
431 rop[55] - rop[56] - rop[57] - rop[58] - rop[59] - rop[60] - rop[
432 61] - rop[62] - rop[63] - rop[64] + rop[66] + rop[72] + rop[73] +
433 rop[74] + rop[76] + rop[80] + rop[84] + rop[98] + rop[100] + rop[
434 101] + rop[102] + rop[104] + rop[107] + rop[109] + rop[110] + rop[
435 112] + rop[117] + rop[119] + rop[121] + rop[131] + rop[136] + rop[
436 145] + rop[146] + rop[150] + rop[150] + rop[152] + rop[152] - rop[
437 156] - rop[160] - rop[161] - rop[163] + rop[164];
438 rrsp[3] = -rop[0] - rop[0] - rop[1] - rop[2] - rop[3] - rop[4] - rop[5] -
439 rop[6] - rop[7] - rop[8] - rop[9] - rop[10] - rop[11] - rop[12] -
440 rop[13] - rop[14] - rop[15] - rop[16] - rop[17] - rop[18] - rop[
441 19] - rop[20] - rop[21] - rop[22] - rop[23] - rop[24] + rop[25] +
442 rop[31] + rop[37] + rop[68] + rop[97] + rop[99] + rop[127] - rop[
443 145] - rop[146] + rop[151] + rop[154] - rop[164] - rop[165];
444 rrsp[4] = rop[0] + rop[3] - rop[25] - rop[26] - rop[27] - rop[31] + rop[
445 38] + rop[69] + rop[90] + rop[93] - rop[97] - rop[99] - rop[109]
446 - rop[117] - rop[118] - rop[127] - rop[128] - rop[138] - rop[139]
447 - rop[140] - rop[142] - rop[143] - rop[150] - rop[151] - rop[154]
448 - rop[155] - rop[158];
449 rrsp[5] = rop[1] + rop[2] + rop[3] + rop[4] + rop[10] + rop[12] + rop[14]
450 + rop[15] + rop[21] + rop[23] + rop[31] - rop[36] + rop[39] + rop[
451 39] + rop[41] + rop[52] - rop[66] - rop[67] - rop[67] - rop[68] -
452 rop[68] - rop[69] - rop[70] - rop[72] - rop[73] - rop[74] - rop[
453 75] - rop[76] - rop[77] - rop[78] - rop[79] - rop[80] - rop[81] -
454 rop[82] - rop[83] - rop[84] - rop[85] - rop[86] - rop[87] - rop[
455 88] - rop[89] + rop[92] + rop[94] + rop[95] + rop[109] + rop[117]
456 + rop[128] + rop[143] - rop[148] + rop[158] - rop[162];
457 rrsp[6] = rop[36] + rop[37] + rop[41] + rop[53] + rop[66] + rop[68] + rop[
458 69] + rop[70] + rop[75] + rop[77] + rop[78] + rop[79] + rop[81] +
459 rop[82] + rop[83] + rop[86] + rop[87] + rop[88] + rop[89] - rop[
460 101] + rop[118] - rop[153] + rop[162];
461 rrsp[7] = -rop[3] + rop[4] + rop[26] + rop[27] - rop[37] - rop[38] - rop[
462 39] + rop[40] - rop[69] + rop[70] - rop[90] - rop[90] - rop[92] -
463 rop[93] - rop[94] - rop[95] - rop[96] + rop[129] + rop[138] + rop[
464 139] + rop[142] + rop[155];
465 rrsp[8] = -rop[4] - rop[40] - rop[41] + rop[67] - rop[70] + rop[90] + rop[
466 96] - rop[129];
467 rrsp[9] = -rop[9] + rop[10] + rop[19] + rop[20] + rop[43] - rop[45] + rop[
468 46] + rop[52] + rop[64] - rop[77] - rop[78] + rop[79] + rop[85] -
469 rop[93] - rop[94] - rop[98] + rop[110] - rop[112] + rop[113] +
470 rop[113] + rop[119] - rop[121] + rop[122] + rop[122] - rop[127] -
471 rop[128] - rop[129] - rop[130] - rop[130] - rop[131] - rop[131] -
472 rop[132] - rop[133] - rop[134] - rop[135] - rop[145] - rop[148] +
473 rop[149] + rop[160] + rop[163] + rop[164] - rop[166];
474 rrsp[10] = -rop[10] + rop[45] - rop[46] - rop[79] + rop[93] - rop[104] -
475 rop[113] - rop[122] + rop[129] + rop[132] + rop[133] + rop[134] +
476 rop[135];
477 rrsp[11] = rop[5] + rop[7] - rop[11] + rop[12] + rop[17] + rop[22] + rop[
478 22] - rop[25] + rop[48] + rop[62] + rop[64] - rop[65] + rop[72] -
479 rop[80] + rop[81] + rop[85] - rop[95] + rop[97] - rop[105] + rop[
480 106] + rop[108] + rop[109] - rop[114] + rop[117] + rop[118] + rop[
481 125] + rop[132] + rop[136] + rop[138] + rop[143] + rop[143] + rop[
482 144] + rop[144] + rop[145] + rop[158];
483 rrsp[12] = rop[11] + rop[13] + rop[24] + rop[25] + rop[80] + rop[95] -
484 rop[106] - rop[125] + rop[150];
485 rrsp[13] = rop[9] - rop[14] + rop[15] + rop[20] - rop[26] + rop[47] - rop[
486 49] - rop[50] + rop[51] + rop[65] + rop[74] + rop[76] - rop[82] +
487 rop[83] + rop[92] - rop[96] + rop[101] - rop[107] + rop[125] +
488 rop[128] - rop[133] + rop[139] + rop[140] + rop[148] + rop[151] +
489 rop[153] + rop[158];
490 rrsp[14] = -rop[16] - rop[17] - rop[54] + rop[56] - rop[84] - rop[85] +
491 rop[86] + rop[98] + rop[102] + rop[108] + rop[111] + rop[141] +
492 rop[144] + rop[152] + rop[155];
493 rrsp[15] = -rop[19] + rop[55] - rop[57] - rop[58] + rop[60] - rop[87] +
494 rop[104] + rop[112] + rop[121] - rop[134] - rop[141] + rop[142] -
495 rop[146] + rop[163];
496 rrsp[16] = -rop[21] + rop[59] - rop[61] - rop[88] + rop[130] - rop[135];
497 rrsp[17] = rop[18] - rop[23] - rop[24] - rop[63] - rop[64] + rop[84] -
498 rop[89] + rop[107] + rop[114] - rop[156] + rop[161] + rop[162] +
499 rop[164];
500 rrsp[18] = -rop[163] - rop[164] - rop[165] + rop[166];
501 rrsp[19] = 0.f;
502
503
504 // convert units from mol/cm3*s to kmol/m3*s
505
506 for(int k=1; k<=19; k++)
507 rrsp[k] *= 1000.0;
508
509}
510
511
512/* ----------------------------------------------------------------------C */
513
514/* Subroutine */ int ratt_(double *t, double *alogt, double *rf,
515 double *rb, double *rklow)
516{
517 /* Initialized data */
518
519 static double small = 1e-50;
520
521 /* Builtin functions */
522 //dol double exp(double);
523
524 /* Local variables */
525 static int n;
526 static double eg[29], ti, ti2, eqk[167], smh[29], tmp, pfac, pfac2,
527 pfac3;
528 extern /* Subroutine */ int rdsmh_(double *, double *, double
529 *);
530
531
532/* *****precision > double */
533/* *****END precision > double */
534
535 /* Parameter adjustments */
536 --rklow;
537 --rb;
538 --rf;
539
540 /* Function Body */
541
542 ti = 1. / *t;
543 ti2 = ti * ti;
544
545 rf[1] = ti * 1.2e17;
546 rf[2] = ti * 5e17;
547 rf[3] = exp(*alogt * 2.7 + 10.5635949 - ti * 3150.13634);
548 rf[4] = 2e13;
549 rf[5] = exp(*alogt * 2. + 16.0803938 - ti * 2012.86667);
550 rf[6] = 5.7e13;
551 rf[7] = 8e13;
552 rf[8] = 1.5e13;
553 rf[9] = 1.5e13;
554 rf[10] = 5.06e13;
555 rf[11] = exp(*alogt * 1.5 + 20.7430685 - ti * 4327.66334);
556 rf[12] = exp(23.6136376 - ti * 1200.17175);
557 rf[13] = 3e13;
558 rf[14] = 3e13;
559 rf[15] = exp(31.2945828 - ti * 1781.387);
560 rf[16] = 1e13;
561 tmp = exp(*alogt * 2. - ti * 956.111669);
562 rf[17] = tmp * 1.35e7;
563 rf[18] = tmp * 6.94e6;
564 rf[19] = 3e13;
565 tmp = exp(*alogt * 1.83 - ti * 110.707667);
566 rf[20] = tmp * 2.5e7;
567 rf[147] = tmp * 3.35e6;
568 rf[21] = 2.24e13;
569 rf[22] = exp(*alogt * 1.92 + 18.3130955 - ti * 2863.30284);
570 rf[23] = 1e14;
571 tmp = exp(ti * -4025.73334);
572 rf[24] = tmp * 1e13;
573 rf[64] = tmp * 5e13;
574 rf[25] = exp(28.1906369 - ti * 679.342501);
575 rf[26] = exp(28.5473118 - ti * 24053.7567);
576 rf[27] = exp(32.2361913 - ti * 20128.6667);
577 rf[28] = exp(42.4761511 - *alogt * .86);
578 tmp = exp(*alogt * -1.24);
579 rf[29] = tmp * 2.08e19;
580 rf[31] = tmp * 2.6e19;
581 rf[30] = exp(43.8677883 - *alogt * .76);
582 rf[32] = exp(37.8159211 - *alogt * .6707 - ti * 8575.31523);
583 rf[33] = ti * 1e18;
584 rf[34] = exp(39.0385861 - *alogt * .6);
585 rf[35] = exp(45.5408762 - *alogt * 1.25);
586 rf[36] = ti2 * 5.5e20;
587 rf[37] = ti2 * 2.2e22;
588 rf[38] = exp(29.0097872 - ti * 337.658384);
589 rf[39] = exp(31.9812991 - ti * 537.435401);
590 rf[40] = exp(31.3686907 - ti * 319.542584);
591 rf[41] = exp(*alogt * 2. + 16.308716 - ti * 2616.72667);
592 rf[42] = exp(29.9336062 - ti * 1811.58);
593 rf[43] = 1.65e14;
594 rf[44] = 6e14;
595 rf[45] = 3e13;
596 rf[46] = exp(37.1706652 - *alogt * .534 - ti * 269.724134);
597 rf[47] = exp(*alogt * 1.62 + 20.3077504 - ti * 5454.86868);
598 rf[48] = exp(*alogt * .48 + 27.7171988 + ti * 130.836334);
599 rf[49] = 7.34e13;
600 rf[50] = exp(*alogt * .454 + 27.014835 - ti * 1308.36334);
601 rf[51] = exp(*alogt * 1.9 + 17.8655549 - ti * 1379.8201);
602 rf[52] = 2e13;
603 rf[53] = exp(*alogt * .5 + 28.0364862 + ti * 55.3538334);
604 rf[54] = exp(33.1993656 - *alogt * .23 - ti * 538.441834);
605 rf[55] = exp(29.3537877 - ti * 1207.72);
606 rf[56] = exp(*alogt * .27 + 29.4360258 - ti * 140.900667);
607 rf[57] = 3.03e13;
608 rf[58] = exp(*alogt * .454 + 27.3863985 - ti * 915.854335);
609 rf[59] = exp(*alogt * 2.53 + 14.096923 - ti * 6159.37201);
610 rf[60] = exp(40.7945264 - *alogt * .99 - ti * 795.082335);
611 rf[61] = 2e12;
612 rf[62] = exp(*alogt * 1.9 + 18.5604427 - ti * 3789.22151);
613 rf[63] = 1e14;
614 rf[65] = exp(30.0558238 - ti * 1725.02674);
615 rf[66] = exp(*alogt * 1.5 + 17.5767107 - ti * 40056.0467);
616 rf[67] = exp(*alogt * 1.51 + 19.190789 - ti * 1726.03317);
617 rf[68] = exp(31.9350862 - *alogt * .37);
618 rf[69] = exp(*alogt * 2.4 + 10.482906 + ti * 1061.78717);
619 rf[70] = exp(ti * 251.608334 + 30.3051698);
620 rf[71] = exp(28.3241683 - ti * 214.873517);
621 rf[72] = exp(41.9771599 - ti * 14799.6022);
622 rf[73] = 5e13;
623 rf[74] = 3e13;
624 rf[75] = 2e13;
625 rf[76] = exp(*alogt * 2. + 16.2403133 - ti * 1509.65);
626 rf[77] = 3e13;
627 rf[78] = exp(*alogt * 1.6 + 17.8408622 - ti * 2727.43434);
628 rf[79] = exp(41.0064751 - *alogt * 1.34 - ti * 713.058018);
629 rf[80] = exp(*alogt * 1.6 + 18.4206807 - ti * 1570.036);
630 rf[81] = exp(*alogt * 1.228 + 17.6783433 - ti * 35.2251667);
631 rf[82] = 5e13;
632 rf[83] = exp(*alogt * 1.18 + 21.9558261 + ti * 224.93785);
633 rf[84] = 5e12;
634 rf[85] = exp(*alogt * 4.5 - 8.4310155 + ti * 503.216668);
635 rf[86] = exp(*alogt * 4. - 7.6354939 + ti * 1006.43334);
636 rf[87] = 5e12;
637 rf[88] = exp(*alogt * 2. + 14.4032972 - ti * 1258.04167);
638 rf[89] = exp(*alogt * 2.12 + 15.0796373 - ti * 437.798501);
639 rf[90] = exp(29.6459241 - ti * 1006.43334);
640 rf[91] = exp(ti * 820.243168 + 25.5908003);
641 rf[92] = exp(33.6712758 - ti * 6038.60001);
642 rf[93] = 2e13;
643 rf[94] = 1e12;
644 rf[95] = 2.87e13;
645 rf[96] = exp(32.6416564 - ti * 11875.9134);
646 rf[97] = exp(*alogt * 2. + 15.5382772 - ti * 6038.60001);
647 rf[98] = exp(31.6914641 - ti * 289.852801);
648 rf[99] = 5e13;
649 rf[100] = 6.71e13;
650 rf[101] = exp(32.3131523 - ti * 1565.00384);
651 rf[102] = exp(ti * 379.928584 + 29.3732401);
652 rf[103] = 4e13;
653 rf[105] = 6e13;
654 rf[106] = 5e13;
655 rf[107] = exp(32.8780452 - ti * 7946.79762);
656 rf[108] = exp(ti * 259.156584 + 32.1806786);
657 rf[109] = 5e13;
658 tmp = exp(ti * -754.825001);
659 rf[110] = tmp * 5e12;
660 rf[151] = tmp * 5.8e12;
661 rf[152] = tmp * 2.4e12;
662 rf[111] = exp(*alogt * 2. + 13.1223634 - ti * 3638.25651);
663 rf[112] = exp(35.00878 - ti * 6010.41988);
664 rf[113] = 4e13;
665 rf[114] = exp(*alogt * 2. + 14.7156719 - ti * 4161.60184);
666 rf[115] = exp(*alogt * .5 + 27.4203001 - ti * 2269.50717);
667 rf[117] = exp(30.3390713 - ti * 301.930001);
668 rf[118] = 2.8e13;
669 rf[119] = 1.2e13;
670 rf[120] = 7e13;
671 rf[121] = 3e13;
672 tmp = exp(ti * 286.833501);
673 rf[122] = tmp * 1.2e13;
674 rf[123] = tmp * 1.6e13;
675 rf[124] = 9e12;
676 rf[125] = 7e12;
677 rf[126] = 1.4e13;
678 rf[128] = exp(31.2033668 - ti * 15338.044);
679 rf[129] = exp(28.4682686 - ti * 10222.8466);
680 rf[130] = exp(*alogt * 2.47 + 10.1064284 - ti * 2606.66234);
681 rf[131] = exp(38.7538626 - *alogt * 1.18 - ti * 329.103701);
682 rf[132] = exp(*alogt * .1 + 29.5538088 - ti * 5334.09668);
683 rf[133] = 2.648e13;
684 rf[134] = exp(*alogt * 2.81 + 8.10772006 - ti * 2948.84967);
685 rf[135] = exp(*alogt * 2. + 12.3327053 - ti * 4629.59334);
686 rf[136] = exp(*alogt * 1.74 + 15.6303353 - ti * 5258.61418);
687 tmp = exp(*alogt * -1. - ti * 8554.68335);
688 rf[137] = tmp * 1.5e18;
689 rf[138] = tmp * 1.87e17;
690 rf[139] = exp(30.2300002 - ti * 201.286667);
691 rf[140] = exp(*alogt * 7.6 - 28.4796532 + ti * 1776.35484);
692 rf[141] = exp(38.3630605 - *alogt * 1.39 - ti * 510.764918);
693 rf[142] = exp(*alogt * .44 + 29.7104627 - ti * 43664.1103);
694 rf[143] = exp(27.4566677 - ti * 1949.96459);
695 rf[144] = exp(28.7941719 - ti * 429.747034);
696 rf[145] = 1e13;
697 rf[146] = 3.37e13;
698 rf[148] = exp(36.1482143 - ti * 8720.74485);
699 rf[149] = exp(*alogt * .5 + 22.8027074 + ti * 883.145252);
700 rf[150] = exp(*alogt * .43 + 28.3090547 + ti * 186.190167);
701 rf[154] = exp(*alogt * .25 + 24.9457104 + ti * 470.507584);
702 rf[155] = exp(*alogt * .29 + 25.5207079 - ti * 5.53538334);
703 rf[156] = exp(*alogt * 1.61 + 14.1059389 + ti * 193.2352);
704 rf[157] = exp(*alogt * .422 + 26.9105027 + ti * 883.145252);
705 rf[159] = 1.81e10;
706 rf[161] = 2.2e13;
707 rf[162] = 1.1e13;
708 rf[163] = 1.2e13;
709 rf[164] = exp(51.6031099 - *alogt * 2.39 - ti * 5625.96234);
710 rf[165] = exp(*alogt * 1.65 + 18.6030023 - ti * 164.55185);
711 rf[166] = exp(*alogt * 1.65 + 17.3708586 + ti * 489.126601);
712 rf[167] = 2.5e13;
713
714 rdsmh_(t, alogt, smh);
715 for (n = 1; n <= 28; ++n) {
716 eg[n - 1] = exp(smh[n - 1]);
717 }
718
719 pfac = 1013250. / (*t * 83145100.);
720 pfac2 = pfac * pfac;
721 pfac3 = pfac2 * pfac;
722
723 eqk[0] = eg[3] / eg[2] / eg[2] / pfac;
724 eqk[1] = eg[4] / eg[1] / eg[2] / pfac;
725 eqk[2] = eg[1] * eg[4] / eg[0] / eg[2];
726 eqk[3] = eg[3] * eg[4] / eg[2] / eg[6];
727 eqk[4] = eg[4] * eg[6] / eg[2] / eg[7];
728 eqk[5] = eg[1] * eg[14] / eg[2] / eg[9];
729 eqk[6] = eg[1] * eg[16] / eg[2] / eg[10];
730 eqk[7] = eg[0] * eg[14] / eg[2] / eg[11];
731 eqk[8] = eg[1] * eg[16] / eg[2] / eg[11];
732 eqk[9] = eg[1] * eg[17] / eg[2] / eg[12];
733 eqk[10] = eg[4] * eg[12] / eg[2] / eg[13];
734 eqk[11] = eg[15] / eg[2] / eg[14] / pfac;
735 eqk[12] = eg[4] * eg[14] / eg[2] / eg[16];
736 eqk[13] = eg[1] * eg[15] / eg[2] / eg[16];
737 eqk[14] = eg[4] * eg[16] / eg[2] / eg[17];
738 eqk[15] = eg[4] * eg[17] / eg[2] / eg[18];
739 eqk[16] = eg[1] * eg[24] / eg[2] / eg[19];
740 eqk[17] = eg[10] * eg[14] / eg[2] / eg[19];
741 eqk[18] = eg[1] * eg[25] / eg[2] / eg[20];
742 eqk[19] = eg[12] * eg[16] / eg[2] / eg[21];
743 eqk[20] = eg[12] * eg[17] / eg[2] / eg[22];
744 eqk[21] = eg[4] * eg[22] / eg[2] / eg[23];
745 eqk[22] = eg[1] * eg[14] * eg[14] / eg[2] / eg[24] * pfac;
746 eqk[23] = eg[4] * eg[24] / eg[2] / eg[25];
747 eqk[24] = eg[10] * eg[15] / eg[2] / eg[25];
748 eqk[25] = eg[2] * eg[15] / eg[3] / eg[14];
749 eqk[26] = eg[6] * eg[16] / eg[3] / eg[17];
750 eqk[27] = eg[6] / eg[1] / eg[3] / pfac;
751 eqk[28] = eqk[27];
752 eqk[29] = eqk[27];
753 eqk[30] = eqk[27];
754 eqk[31] = eg[2] * eg[4] / eg[1] / eg[3];
755 eqk[32] = eg[0] / eg[1] / eg[1] / pfac;
756 eqk[33] = eqk[32];
757 eqk[34] = eqk[32];
758 eqk[35] = eqk[32];
759 eqk[36] = eg[5] / eg[1] / eg[4] / pfac;
760 eqk[37] = eg[2] * eg[5] / eg[1] / eg[6];
761 eqk[38] = eg[0] * eg[3] / eg[1] / eg[6];
762 eqk[39] = eg[4] * eg[4] / eg[1] / eg[6];
763 eqk[40] = eg[0] * eg[6] / eg[1] / eg[7];
764 eqk[41] = eg[4] * eg[5] / eg[1] / eg[7];
765 eqk[42] = eg[0] * eg[8] / eg[1] / eg[9];
766 eqk[43] = eg[12] / eg[1] / eg[10] / pfac;
767 eqk[44] = eg[0] * eg[9] / eg[1] / eg[11];
768 eqk[45] = eg[13] / eg[1] / eg[12] / pfac;
769 eqk[46] = eg[0] * eg[12] / eg[1] / eg[13];
770 eqk[47] = eg[17] / eg[1] / eg[16] / pfac;
771 eqk[48] = eg[0] * eg[14] / eg[1] / eg[16];
772 eqk[49] = eg[18] / eg[1] / eg[17] / pfac;
773 eqk[50] = eg[0] * eg[16] / eg[1] / eg[17];
774 eqk[51] = eg[0] * eg[17] / eg[1] / eg[18];
775 eqk[52] = eg[4] * eg[12] / eg[1] / eg[18];
776 eqk[53] = eg[5] * eg[11] / eg[1] / eg[18];
777 eqk[54] = eg[20] / eg[1] / eg[19] / pfac;
778 eqk[55] = eg[21] / eg[1] / eg[20] / pfac;
779 eqk[56] = eg[0] * eg[19] / eg[1] / eg[20];
780 eqk[57] = eg[22] / eg[1] / eg[21] / pfac;
781 eqk[58] = eg[0] * eg[20] / eg[1] / eg[21];
782 eqk[59] = eg[23] / eg[1] / eg[22] / pfac;
783 eqk[60] = eg[0] * eg[21] / eg[1] / eg[22];
784 eqk[61] = eg[0] * eg[22] / eg[1] / eg[23];
785 eqk[62] = eg[11] * eg[14] / eg[1] / eg[24];
786 eqk[63] = eg[0] * eg[24] / eg[1] / eg[25];
787 eqk[64] = eg[12] * eg[14] / eg[1] / eg[25];
788 eqk[65] = eg[17] / eg[0] / eg[14] / pfac;
789 eqk[66] = eg[1] * eg[5] / eg[0] / eg[4];
790 eqk[67] = eg[7] / eg[4] / eg[4] / pfac;
791 eqk[68] = eg[2] * eg[5] / eg[4] / eg[4];
792 eqk[69] = eg[3] * eg[5] / eg[4] / eg[6];
793 eqk[147] = eqk[69];
794 eqk[70] = eg[5] * eg[6] / eg[4] / eg[7];
795 eqk[71] = eqk[70];
796 eqk[72] = eg[1] * eg[14] / eg[4] / eg[8];
797 eqk[73] = eg[1] * eg[16] / eg[4] / eg[9];
798 eqk[74] = eg[1] * eg[17] / eg[4] / eg[10];
799 eqk[75] = eg[5] * eg[9] / eg[4] / eg[10];
800 eqk[76] = eg[1] * eg[17] / eg[4] / eg[11];
801 eqk[77] = eg[5] * eg[10] / eg[4] / eg[12];
802 eqk[78] = eg[5] * eg[11] / eg[4] / eg[12];
803 eqk[79] = eg[5] * eg[12] / eg[4] / eg[13];
804 eqk[80] = eg[1] * eg[15] / eg[4] / eg[14];
805 eqk[81] = eg[5] * eg[14] / eg[4] / eg[16];
806 eqk[82] = eg[5] * eg[16] / eg[4] / eg[17];
807 eqk[83] = eg[5] * eg[17] / eg[4] / eg[18];
808 eqk[84] = eg[1] * eg[25] / eg[4] / eg[19];
809 eqk[85] = eg[12] * eg[14] / eg[4] / eg[19];
810 eqk[86] = eg[5] * eg[19] / eg[4] / eg[20];
811 eqk[87] = eg[5] * eg[20] / eg[4] / eg[21];
812 eqk[88] = eg[5] * eg[22] / eg[4] / eg[23];
813 eqk[89] = eg[5] * eg[24] / eg[4] / eg[25];
814 eqk[90] = eg[3] * eg[7] / eg[6] / eg[6];
815 eqk[91] = eqk[90];
816 eqk[92] = eg[4] * eg[17] / eg[6] / eg[10];
817 eqk[93] = eg[3] * eg[13] / eg[6] / eg[12];
818 eqk[94] = eg[4] * eg[18] / eg[6] / eg[12];
819 eqk[95] = eg[4] * eg[15] / eg[6] / eg[14];
820 eqk[96] = eg[7] * eg[16] / eg[6] / eg[17];
821 eqk[97] = eg[2] * eg[14] / eg[3] / eg[8];
822 eqk[98] = eg[1] * eg[19] / eg[8] / eg[12];
823 eqk[99] = eg[2] * eg[16] / eg[3] / eg[9];
824 eqk[100] = eg[1] * eg[10] / eg[0] / eg[9];
825 eqk[101] = eg[1] * eg[17] / eg[5] / eg[9];
826 eqk[102] = eg[1] * eg[19] / eg[9] / eg[10];
827 eqk[104] = eg[1] * eg[21] / eg[9] / eg[13];
828 eqk[105] = eg[24] / eg[9] / eg[14] / pfac;
829 eqk[106] = eg[14] * eg[16] / eg[9] / eg[15];
830 eqk[107] = eg[1] * eg[25] / eg[9] / eg[17];
831 eqk[108] = eg[14] * eg[19] / eg[9] / eg[24];
832 eqk[110] = eg[1] * eg[12] / eg[0] / eg[10];
833 eqk[111] = eg[0] * eg[19] / eg[10] / eg[10];
834 eqk[112] = eg[1] * eg[21] / eg[10] / eg[12];
835 eqk[113] = eg[12] * eg[12] / eg[10] / eg[13];
836 eqk[114] = eg[25] / eg[10] / eg[14] / pfac;
837 eqk[116] = eg[10] / eg[11];
838 eqk[120] = eqk[116];
839 eqk[123] = eqk[116];
840 eqk[124] = eqk[116];
841 eqk[117] = eg[1] * eg[4] * eg[14] / eg[3] / eg[11] * pfac;
842 eqk[118] = eg[5] * eg[14] / eg[3] / eg[11];
843 eqk[119] = eg[1] * eg[12] / eg[0] / eg[11];
844 eqk[121] = eg[1] * eg[21] / eg[11] / eg[12];
845 eqk[122] = eg[12] * eg[12] / eg[11] / eg[13];
846 eqk[125] = eg[14] * eg[17] / eg[11] / eg[15];
847 eqk[127] = eg[2] * eg[18] / eg[3] / eg[12];
848 eqk[128] = eg[4] * eg[17] / eg[3] / eg[12];
849 eqk[129] = eg[6] * eg[13] / eg[7] / eg[12];
850 eqk[130] = eg[23] / eg[12] / eg[12] / pfac;
851 eqk[131] = eg[1] * eg[22] / eg[12] / eg[12];
852 eqk[132] = eg[13] * eg[14] / eg[12] / eg[16];
853 eqk[133] = eg[13] * eg[16] / eg[12] / eg[17];
854 eqk[134] = eg[13] * eg[20] / eg[12] / eg[21];
855 eqk[135] = eg[13] * eg[22] / eg[12] / eg[23];
856 eqk[136] = eg[1] * eg[14] / eg[16] * pfac;
857 eqk[137] = eqk[136];
858 eqk[138] = eg[6] * eg[14] / eg[3] / eg[16];
859 eqk[139] = eg[6] * eg[17] / eg[3] / eg[18];
860 eqk[140] = eg[16] * eg[17] / eg[3] / eg[20];
861 eqk[141] = eg[0] * eg[19] / eg[21] * pfac;
862 eqk[142] = eg[6] * eg[21] / eg[3] / eg[22];
863 eqk[143] = eg[4] * eg[14] * eg[14] / eg[3] / eg[24] * pfac;
864 eqk[144] = eg[14] * eg[14] * eg[19] / eg[24] / eg[24] * pfac;
865 eqk[146] = eg[1] * eg[26] / eg[2] / eg[21];
866 eqk[149] = eg[12] / eg[0] / eg[9] / pfac;
867 eqk[151] = eg[2] * eg[17] / eg[3] / eg[10];
868 eqk[155] = eg[6] * eg[19] / eg[3] / eg[20];
869 eqk[156] = eg[26] / eg[1] / eg[25] / pfac;
870 eqk[160] = eg[12] * eg[16] / eg[1] / eg[26];
871 eqk[161] = eg[0] * eg[25] / eg[1] / eg[26];
872 eqk[162] = eg[5] * eg[25] / eg[4] / eg[26];
873 eqk[163] = eg[12] * eg[21] / eg[1] / eg[27];
874 eqk[164] = eg[1] * eg[12] * eg[25] / eg[2] / eg[27] * pfac;
875 eqk[166] = eg[27] / eg[12] / eg[20] / pfac;
876
877 rb[1] = rf[1] / max(eqk[0],small);
878 rb[2] = rf[2] / max(eqk[1],small);
879 rb[3] = rf[3] / max(eqk[2],small);
880 rb[4] = rf[4] / max(eqk[3],small);
881 rb[5] = rf[5] / max(eqk[4],small);
882 rb[6] = rf[6] / max(eqk[5],small);
883 rb[7] = rf[7] / max(eqk[6],small);
884 rb[8] = rf[8] / max(eqk[7],small);
885 rb[9] = rf[9] / max(eqk[8],small);
886 rb[10] = rf[10] / max(eqk[9],small);
887 rb[11] = rf[11] / max(eqk[10],small);
888 rb[12] = rf[12] / max(eqk[11],small);
889 rb[13] = rf[13] / max(eqk[12],small);
890 rb[14] = rf[14] / max(eqk[13],small);
891 rb[15] = rf[15] / max(eqk[14],small);
892 rb[16] = rf[16] / max(eqk[15],small);
893 rb[17] = rf[17] / max(eqk[16],small);
894 rb[18] = rf[18] / max(eqk[17],small);
895 rb[19] = rf[19] / max(eqk[18],small);
896 rb[20] = rf[20] / max(eqk[19],small);
897 rb[21] = rf[21] / max(eqk[20],small);
898 rb[22] = rf[22] / max(eqk[21],small);
899 rb[23] = rf[23] / max(eqk[22],small);
900 rb[24] = rf[24] / max(eqk[23],small);
901 rb[25] = rf[25] / max(eqk[24],small);
902 rb[26] = rf[26] / max(eqk[25],small);
903 rb[27] = rf[27] / max(eqk[26],small);
904 rb[28] = rf[28] / max(eqk[27],small);
905 rb[29] = rf[29] / max(eqk[28],small);
906 rb[30] = rf[30] / max(eqk[29],small);
907 rb[31] = rf[31] / max(eqk[30],small);
908 rb[32] = rf[32] / max(eqk[31],small);
909 rb[33] = rf[33] / max(eqk[32],small);
910 rb[34] = rf[34] / max(eqk[33],small);
911 rb[35] = rf[35] / max(eqk[34],small);
912 rb[36] = rf[36] / max(eqk[35],small);
913 rb[37] = rf[37] / max(eqk[36],small);
914 rb[38] = rf[38] / max(eqk[37],small);
915 rb[39] = rf[39] / max(eqk[38],small);
916 rb[40] = rf[40] / max(eqk[39],small);
917 rb[41] = rf[41] / max(eqk[40],small);
918 rb[42] = rf[42] / max(eqk[41],small);
919 rb[43] = rf[43] / max(eqk[42],small);
920 rb[44] = rf[44] / max(eqk[43],small);
921 rb[45] = rf[45] / max(eqk[44],small);
922 rb[46] = rf[46] / max(eqk[45],small);
923 rb[47] = rf[47] / max(eqk[46],small);
924 rb[48] = rf[48] / max(eqk[47],small);
925 rb[49] = rf[49] / max(eqk[48],small);
926 rb[50] = rf[50] / max(eqk[49],small);
927 rb[51] = rf[51] / max(eqk[50],small);
928 rb[52] = rf[52] / max(eqk[51],small);
929 rb[53] = rf[53] / max(eqk[52],small);
930 rb[54] = rf[54] / max(eqk[53],small);
931 rb[55] = rf[55] / max(eqk[54],small);
932 rb[56] = rf[56] / max(eqk[55],small);
933 rb[57] = rf[57] / max(eqk[56],small);
934 rb[58] = rf[58] / max(eqk[57],small);
935 rb[59] = rf[59] / max(eqk[58],small);
936 rb[60] = rf[60] / max(eqk[59],small);
937 rb[61] = rf[61] / max(eqk[60],small);
938 rb[62] = rf[62] / max(eqk[61],small);
939 rb[63] = rf[63] / max(eqk[62],small);
940 rb[64] = rf[64] / max(eqk[63],small);
941 rb[65] = rf[65] / max(eqk[64],small);
942 rb[66] = rf[66] / max(eqk[65],small);
943 rb[67] = rf[67] / max(eqk[66],small);
944 rb[68] = rf[68] / max(eqk[67],small);
945 rb[69] = rf[69] / max(eqk[68],small);
946 rb[70] = rf[70] / max(eqk[69],small);
947 rb[71] = rf[71] / max(eqk[70],small);
948 rb[72] = rf[72] / max(eqk[71],small);
949 rb[73] = rf[73] / max(eqk[72],small);
950 rb[74] = rf[74] / max(eqk[73],small);
951 rb[75] = rf[75] / max(eqk[74],small);
952 rb[76] = rf[76] / max(eqk[75],small);
953 rb[77] = rf[77] / max(eqk[76],small);
954 rb[78] = rf[78] / max(eqk[77],small);
955 rb[79] = rf[79] / max(eqk[78],small);
956 rb[80] = rf[80] / max(eqk[79],small);
957 rb[81] = rf[81] / max(eqk[80],small);
958 rb[82] = rf[82] / max(eqk[81],small);
959 rb[83] = rf[83] / max(eqk[82],small);
960 rb[84] = rf[84] / max(eqk[83],small);
961 rb[85] = rf[85] / max(eqk[84],small);
962 rb[86] = rf[86] / max(eqk[85],small);
963 rb[87] = rf[87] / max(eqk[86],small);
964 rb[88] = rf[88] / max(eqk[87],small);
965 rb[89] = rf[89] / max(eqk[88],small);
966 rb[90] = rf[90] / max(eqk[89],small);
967 rb[91] = rf[91] / max(eqk[90],small);
968 rb[92] = rf[92] / max(eqk[91],small);
969 rb[93] = rf[93] / max(eqk[92],small);
970 rb[94] = rf[94] / max(eqk[93],small);
971 rb[95] = rf[95] / max(eqk[94],small);
972 rb[96] = rf[96] / max(eqk[95],small);
973 rb[97] = rf[97] / max(eqk[96],small);
974 rb[98] = rf[98] / max(eqk[97],small);
975 rb[99] = rf[99] / max(eqk[98],small);
976 rb[100] = rf[100] / max(eqk[99],small);
977 rb[101] = rf[101] / max(eqk[100],small);
978 rb[102] = rf[102] / max(eqk[101],small);
979 rb[103] = rf[103] / max(eqk[102],small);
980 rf[103] = 0.;
981 rb[105] = rf[105] / max(eqk[104],small);
982 rb[106] = rf[106] / max(eqk[105],small);
983 rb[107] = rf[107] / max(eqk[106],small);
984 rb[108] = rf[108] / max(eqk[107],small);
985 rb[109] = rf[109] / max(eqk[108],small);
986 rf[109] = 0.;
987 rb[110] = 0.;
988 rb[111] = rf[111] / max(eqk[110],small);
989 rb[112] = rf[112] / max(eqk[111],small);
990 rf[112] = 0.;
991 rb[113] = rf[113] / max(eqk[112],small);
992 rb[114] = rf[114] / max(eqk[113],small);
993 rb[115] = rf[115] / max(eqk[114],small);
994 rb[117] = rf[117] / max(eqk[116],small);
995 rb[118] = rf[118] / max(eqk[117],small);
996 rb[119] = rf[119] / max(eqk[118],small);
997 rb[120] = rf[120] / max(eqk[119],small);
998 rb[121] = rf[121] / max(eqk[120],small);
999 rb[122] = rf[122] / max(eqk[121],small);
1000 rb[123] = rf[123] / max(eqk[122],small);
1001 rb[124] = rf[124] / max(eqk[123],small);
1002 rb[125] = rf[125] / max(eqk[124],small);
1003 rb[126] = rf[126] / max(eqk[125],small);
1004 rb[128] = rf[128] / max(eqk[127],small);
1005 rb[129] = rf[129] / max(eqk[128],small);
1006 rb[130] = rf[130] / max(eqk[129],small);
1007 rb[131] = rf[131] / max(eqk[130],small);
1008 rb[132] = rf[132] / max(eqk[131],small);
1009 rb[133] = rf[133] / max(eqk[132],small);
1010 rb[134] = rf[134] / max(eqk[133],small);
1011 rb[135] = rf[135] / max(eqk[134],small);
1012 rb[136] = rf[136] / max(eqk[135],small);
1013 rb[137] = rf[137] / max(eqk[136],small);
1014 rb[138] = rf[138] / max(eqk[137],small);
1015 rb[139] = rf[139] / max(eqk[138],small);
1016 rb[140] = rf[140] / max(eqk[139],small);
1017 rb[142] = rf[142] / max(eqk[141],small);
1018 rb[143] = rf[143] / max(eqk[142],small);
1019 rb[144] = rf[144] / max(eqk[143],small);
1020 rb[145] = rf[145] / max(eqk[144],small);
1021 rf[145] = 0.;
1022 rb[146] = 0.f;
1023 rb[147] = rf[147] / max(eqk[146],small);
1024 rb[148] = rf[148] / max(eqk[147],small);
1025 rb[149] = 0.f;
1026 rb[150] = rf[150] / max(eqk[149],small);
1027 rb[151] = 0.f;
1028 rb[152] = rf[152] / max(eqk[151],small);
1029 rb[153] = 0.f;
1030 rb[154] = 0.f;
1031 rb[156] = rf[156] / max(eqk[155],small);
1032 rb[157] = rf[157] / max(eqk[156],small);
1033 rb[159] = 0.f;
1034 rb[161] = rf[161] / max(eqk[160],small);
1035 rf[161] = 0.;
1036 rb[162] = rf[162] / max(eqk[161],small);
1037 rb[163] = rf[163] / max(eqk[162],small);
1038 rb[164] = rf[164] / max(eqk[163],small);
1039 rb[165] = rf[165] / max(eqk[164],small);
1040 rb[167] = rf[167] / max(eqk[166],small);
1041
1042 rklow[1] = exp(34.0312786 - 1509.65 / *t);
1043 rklow[2] = exp(59.9064331 - *alogt * 2.76 - 805.146668 / *t);
1044 rklow[3] = exp(76.9484824 - *alogt * 4.76 - 1227.84867 / *t);
1045 rklow[4] = exp(56.1662604 - *alogt * 2.57 - 213.867084 / *t);
1046 rklow[5] = exp(69.8660102 - *alogt * 4.8 - 2797.88467 / *t);
1047 rklow[6] = exp(93.4384048 - *alogt * 7.27 - 3633.22434 / *t);
1048 rklow[7] = exp(69.414025 - *alogt * 3.86 - 1670.67934 / *t);
1049 rklow[8] = exp(96.1977483 - *alogt * 7.62 - 3507.42017 / *t);
1050 rklow[9] = exp(95.0941235 - *alogt * 7.08 - 3364.00342 / *t);
1051 rklow[10] = exp(63.7931383 - *alogt * 3.42 - 42446.3259 / *t);
1052 rklow[11] = exp(42.2794408 - *alogt * .9 + 855.468335 / *t);
1053 rklow[12] = exp(65.4619238 - *alogt * 3.74 - 974.227469 / *t);
1054 rklow[13] = exp(76.9748493 - *alogt * 5.11 - 3570.32226 / *t);
1055 rklow[14] = exp(95.6297642 - *alogt * 7.03 - 1389.88444 / *t);
1056 rklow[15] = exp(117.889265 - *alogt * 9.3 - 49214.5901 / *t);
1057 rklow[16] = exp(59.1374013 - *alogt * 2.8 - 296.897834 / *t);
1058 rklow[17] = exp(96.7205025 - *alogt * 7.63 - 1939.39704 / *t);
1059 rklow[18] = exp(135.001549 - *alogt * 11.94 - 4916.3262 / *t);
1060
1061 return 0;
1062} /* ratt_ */
1063
1064/* C */
1065/* ----------------------------------------------------------------------C */
1066/* C */
1067/* Subroutine */ int rdsmh_(double *t, double *tlog, double *smh)
1068{
1069 static double ti, tn[5];
1070
1071
1072/* START PROLOGUE */
1073
1074/* SUBROUTINE DKSMH (T, ICKWRK, RCKWRK, SMH)* */
1075/* Returns the array of entropies minus enthalpies for the species. */
1076/* It is normally not called directly by the user. */
1077
1078/* INPUT */
1079/* T - Temperature. */
1080/* cgs units - K */
1081/* Data type - real scalar */
1082/* TLOG - Log Temperature. */
1083/* OUTPUT */
1084/* SMH - Entropy minus enthalpy for the species, */
1085/* SMH(K) = S(K)/R - H(K)/RT. */
1086/* cgs units - none */
1087/* Data type - real array */
1088/* Dimension SMH(*) at least KK, the total number of */
1089/* species. */
1090
1091/* END PROLOGUE */
1092
1093/* *****precision > double */
1094/* *****END precision > double */
1095/* *****precision > single */
1096/* IMPLICIT REAL (A-H, O-Z), int (I-N) */
1097/* *****END precision > single */
1098
1099
1100 /* Parameter adjustments */
1101 --smh;
1102
1103 /* Function Body */
1104 ti = 1. / *t;
1105
1106 tn[0] = *tlog - 1.f;
1107 tn[1] = *t;
1108 tn[2] = tn[1] * *t;
1109 tn[3] = tn[2] * *t;
1110 tn[4] = tn[3] * *t;
1111
1112 if (*t > 1e3) {
1113
1114 smh[1] = ti * 950.158922 - 3.20502331 + tn[0] * 3.3372792 - tn[1] *
1115 2.47012365e-5 + tn[2] * 8.32427963e-8 - tn[3] *
1116 1.49638662e-11 + tn[4] * 1.00127688e-15;
1117 smh[2] = -.446682914 - ti * 25473.6599 + tn[0] * 2.50000001 - tn[1] *
1118 1.15421486e-11 + tn[2] * 2.69269913e-15 - tn[3] *
1119 3.94596029e-19 + tn[4] * 2.49098679e-23;
1120 smh[3] = 4.78433864 - ti * 29217.5791 + tn[0] * 2.56942078 - tn[1] *
1121 4.29870569e-5 + tn[2] * 6.99140982e-9 - tn[3] *
1122 8.34814992e-13 + tn[4] * 6.14168455e-17;
1123 smh[4] = ti * 1088.45772 + 5.45323129 + tn[0] * 3.28253784 + tn[1] *
1124 7.4154377e-4 - tn[2] * 1.26327778e-7 + tn[3] * 1.74558796e-11
1125 - tn[4] * 1.08358897e-15;
1126 smh[5] = 4.4766961 - ti * 3858.657 + tn[0] * 3.09288767 + tn[1] *
1127 2.74214858e-4 + tn[2] * 2.10842047e-8 - tn[3] * 7.3288463e-12
1128 + tn[4] * 5.8706188e-16;
1129 smh[6] = ti * 30004.2971 + 4.9667701 + tn[0] * 3.03399249 + tn[1] *
1130 .00108845902 - tn[2] * 2.73454197e-8 - tn[3] * 8.08683225e-12
1131 + tn[4] * 8.4100496e-16;
1132 smh[7] = 3.78510215 - ti * 111.856713 + tn[0] * 4.0172109 + tn[1] *
1133 .00111991007 - tn[2] * 1.05609692e-7 + tn[3] * 9.52053083e-12
1134 - tn[4] * 5.39542675e-16;
1135 smh[8] = ti * 17861.7877 + 2.91615662 + tn[0] * 4.16500285 + tn[1] *
1136 .00245415847 - tn[2] * 3.16898708e-7 + tn[3] * 3.09321655e-11
1137 - tn[4] * 1.43954153e-15;
1138 smh[9] = 4.80150373 - ti * 85451.2953 + tn[0] * 2.49266888 + tn[1] *
1139 2.39944642e-5 - tn[2] * 1.20722503e-8 + tn[3] *
1140 3.11909191e-12 - tn[4] * 2.43638946e-16;
1141 smh[10] = 5.48497999 - ti * 71012.4364 + tn[0] * 2.87846473 + tn[1] *
1142 4.85456841e-4 + tn[2] * 2.40742758e-8 - tn[3] *
1143 1.08906541e-11 + tn[4] * 8.80396915e-16;
1144 smh[11] = 6.17119324 - ti * 46263.604 + tn[0] * 2.87410113 + tn[1] *
1145 .00182819646 - tn[2] * 2.34824328e-7 + tn[3] * 2.16816291e-11
1146 - tn[4] * 9.38637835e-16;
1147 smh[12] = 8.62650169 - ti * 50925.9997 + tn[0] * 2.29203842 + tn[1] *
1148 .00232794319 - tn[2] * 3.35319912e-7 + tn[3] * 3.48255e-11 -
1149 tn[4] * 1.69858183e-15;
1150 smh[13] = 8.48007179 - ti * 16775.5843 + tn[0] * 2.28571772 + tn[1] *
1151 .00361995018 - tn[2] * 4.97857247e-7 + tn[3] * 4.9640387e-11
1152 - tn[4] * 2.33577197e-15;
1153 smh[14] = ti * 9468.34459 + 18.437318 + tn[0] * .074851495 + tn[1] *
1154 .00669547335 - tn[2] * 9.55476348e-7 + tn[3] * 1.01910446e-10
1155 - tn[4] * 5.0907615e-15;
1156 smh[15] = ti * 14151.8724 + 7.81868772 + tn[0] * 2.71518561 + tn[1] *
1157 .00103126372 - tn[2] * 1.66470962e-7 + tn[3] * 1.9171084e-11
1158 - tn[4] * 1.01823858e-15;
1159 smh[16] = ti * 48759.166 + 2.27163806 + tn[0] * 3.85746029 + tn[1] *
1160 .00220718513 - tn[2] * 3.69135673e-7 + tn[3] * 4.36241823e-11
1161 - tn[4] * 2.36042082e-15;
1162 smh[17] = 9.79834492 - ti * 4011.91815 + tn[0] * 2.77217438 + tn[1] *
1163 .00247847763 - tn[2] * 4.14076022e-7 + tn[3] * 4.90968148e-11
1164 - tn[4] * 2.66754356e-15;
1165 smh[18] = ti * 13995.8323 + 13.656323 + tn[0] * 1.76069008 + tn[1] *
1166 .00460000041 - tn[2] * 7.37098022e-7 + tn[3] * 8.38676767e-11
1167 - tn[4] * 4.4192782e-15;
1168 smh[19] = 2.929575 - ti * 127.83252 + tn[0] * 3.770799 + tn[1] *
1169 .0039357485 - tn[2] * 4.42730667e-7 + tn[3] * 3.28702583e-11
1170 - tn[4] * 1.056308e-15;
1171 smh[20] = -1.23028121 - ti * 25935.9992 + tn[0] * 4.14756964 + tn[1] *
1172 .00298083332 - tn[2] * 3.9549142e-7 + tn[3] * 3.89510143e-11
1173 - tn[4] * 1.80617607e-15;
1174 smh[21] = 7.78732378 - ti * 34612.8739 + tn[0] * 3.016724 + tn[1] *
1175 .0051651146 - tn[2] * 7.80137248e-7 + tn[3] * 8.480274e-11 -
1176 tn[4] * 4.3130352e-15;
1177 smh[22] = 10.3053693 - ti * 4939.88614 + tn[0] * 2.03611116 + tn[1] *
1178 .00732270755 - tn[2] * 1.11846319e-6 + tn[3] * 1.22685769e-10
1179 - tn[4] * 6.28530305e-15;
1180 smh[23] = 13.4624343 - ti * 12857.52 + tn[0] * 1.95465642 + tn[1] *
1181 .0086986361 - tn[2] * 1.33034445e-6 + tn[3] * 1.46014741e-10
1182 - tn[4] * 7.4820788e-15;
1183 smh[24] = ti * 11426.3932 + 15.1156107 + tn[0] * 1.0718815 + tn[1] *
1184 .0108426339 - tn[2] * 1.67093445e-6 + tn[3] * 1.84510001e-10
1185 - tn[4] * 9.5001445e-15;
1186 smh[25] = -3.9302595 - ti * 19327.215 + tn[0] * 5.6282058 + tn[1] *
1187 .00204267005 - tn[2] * 2.65575783e-7 + tn[3] * 2.38550433e-11
1188 - tn[4] * 9.703916e-16;
1189 smh[26] = ti * 7551.05311 + .632247205 + tn[0] * 4.51129732 + tn[1] *
1190 .00450179872 - tn[2] * 6.94899392e-7 + tn[3] * 7.69454902e-11
1191 - tn[4] * 3.974191e-15;
1192 smh[27] = -5.0320879 - ti * 490.32178 + tn[0] * 5.9756699 + tn[1] *
1193 .0040652957 - tn[2] * 4.5727075e-7 + tn[3] * 3.39192008e-11 -
1194 tn[4] * 1.08800855e-15;
1195 smh[28] = ti * 923.5703 - 13.31335 + tn[0] * 6.732257 + tn[1] *
1196 .00745417 - tn[2] * 8.24983167e-7 + tn[3] * 6.01001833e-11 -
1197 tn[4] * 1.883102e-15;
1198
1199 } else {
1200
1201 smh[1] = ti * 917.935173 + .683010238 + tn[0] * 2.34433112 + tn[1] *
1202 .00399026037 - tn[2] * 3.2463585e-6 + tn[3] * 1.67976745e-9 -
1203 tn[4] * 3.68805881e-13;
1204 smh[2] = -.446682853 - ti * 25473.6599 + tn[0] * 2.5 + tn[1] *
1205 3.52666409e-13 - tn[2] * 3.32653273e-16 + tn[3] *
1206 1.91734693e-19 - tn[4] * 4.63866166e-23;
1207 smh[3] = 2.05193346 - ti * 29122.2592 + tn[0] * 3.1682671 - tn[1] *
1208 .00163965942 + tn[2] * 1.10717733e-6 - tn[3] * 5.10672187e-10
1209 + tn[4] * 1.05632986e-13;
1210 smh[4] = ti * 1063.94356 + 3.65767573 + tn[0] * 3.78245636 - tn[1] *
1211 .00149836708 + tn[2] * 1.641217e-6 - tn[3] * 8.06774591e-10 +
1212 tn[4] * 1.62186419e-13;
1213 smh[5] = -.103925458 - ti * 3615.08056 + tn[0] * 3.99201543 - tn[1] *
1214 .00120065876 + tn[2] * 7.69656402e-7 - tn[3] * 3.23427778e-10
1215 + tn[4] * 6.8205735e-14;
1216 smh[6] = ti * 30293.7267 - .849032208 + tn[0] * 4.19864056 - tn[1] *
1217 .00101821705 + tn[2] * 1.08673369e-6 - tn[3] * 4.57330885e-10
1218 + tn[4] * 8.85989085e-14;
1219 smh[7] = 3.71666245 - ti * 294.80804 + tn[0] * 4.30179801 - tn[1] *
1220 .00237456025 + tn[2] * 3.52638152e-6 - tn[3] * 2.02303245e-9
1221 + tn[4] * 4.64612562e-13;
1222 smh[8] = ti * 17702.5821 + 3.43505074 + tn[0] * 4.27611269 - tn[1] *
1223 2.71411208e-4 + tn[2] * 2.78892835e-6 - tn[3] * 1.79809011e-9
1224 + tn[4] * 4.31227182e-13;
1225 smh[9] = 4.53130848 - ti * 85443.8832 + tn[0] * 2.55423955 - tn[1] *
1226 1.60768862e-4 + tn[2] * 1.22298708e-7 - tn[3] *
1227 6.10195741e-11 + tn[4] * 1.33260723e-14;
1228 smh[10] = 2.08401108 - ti * 70797.2934 + tn[0] * 3.48981665 + tn[1] *
1229 1.61917771e-4 - tn[2] * 2.81498442e-7 + tn[3] *
1230 2.63514439e-10 - tn[4] * 7.03045335e-14;
1231 smh[11] = 1.56253185 - ti * 46004.0401 + tn[0] * 3.76267867 + tn[1] *
1232 4.84436072e-4 + tn[2] * 4.65816402e-7 - tn[3] *
1233 3.20909294e-10 + tn[4] * 8.43708595e-14;
1234 smh[12] = -.769118967 - ti * 50496.8163 + tn[0] * 4.19860411 - tn[1] *
1235 .0011833071 + tn[2] * 1.37216037e-6 - tn[3] * 5.57346651e-10
1236 + tn[4] * 9.71573685e-14;
1237 smh[13] = 1.60456433 - ti * 16444.9988 + tn[0] * 3.6735904 + tn[1] *
1238 .00100547588 + tn[2] * 9.55036427e-7 - tn[3] * 5.72597854e-10
1239 + tn[4] * 1.27192867e-13;
1240 smh[14] = ti * 10246.6476 - 4.64130376 + tn[0] * 5.14987613 - tn[1] *
1241 .0068354894 + tn[2] * 8.19667665e-6 - tn[3] * 4.03952522e-9 +
1242 tn[4] * 8.3346978e-13;
1243 smh[15] = ti * 14344.086 + 3.50840928 + tn[0] * 3.57953347 - tn[1] *
1244 3.0517684e-4 + tn[2] * 1.69469055e-7 + tn[3] * 7.55838237e-11
1245 - tn[4] * 4.52212249e-14;
1246 smh[16] = ti * 48371.9697 + 9.90105222 + tn[0] * 2.35677352 + tn[1] *
1247 .00449229839 - tn[2] * 1.18726045e-6 + tn[3] * 2.04932518e-10
1248 - tn[4] * 7.1849774e-15;
1249 smh[17] = 3.39437243 - ti * 3839.56496 + tn[0] * 4.22118584 - tn[1] *
1250 .00162196266 + tn[2] * 2.29665743e-6 - tn[3] * 1.10953411e-9
1251 + tn[4] * 2.16884433e-13;
1252 smh[18] = ti * 14308.9567 + .6028129 + tn[0] * 4.79372315 - tn[1] *
1253 .00495416685 + tn[2] * 6.22033347e-6 - tn[3] * 3.16071051e-9
1254 + tn[4] * 6.5886326e-13;
1255 smh[19] = 13.152177 - ti * 978.6011 + tn[0] * 2.106204 + tn[1] *
1256 .0036082975 + tn[2] * 8.89745333e-7 - tn[3] * 6.14803e-10 +
1257 tn[4] * 1.037805e-13;
1258 smh[20] = 13.9397051 - ti * 26428.9807 + tn[0] * .808681094 + tn[1] *
1259 .0116807815 - tn[2] * 5.91953025e-6 + tn[3] * 2.33460364e-9 -
1260 tn[4] * 4.25036487e-13;
1261 smh[21] = 8.51054025 - ti * 34859.8468 + tn[0] * 3.21246645 + tn[1] *
1262 7.5739581e-4 + tn[2] * 4.32015687e-6 - tn[3] * 2.98048206e-9
1263 + tn[4] * 7.35754365e-13;
1264 smh[22] = 4.09733096 - ti * 5089.77593 + tn[0] * 3.95920148 - tn[1] *
1265 .00378526124 + tn[2] * 9.51650487e-6 - tn[3] * 5.76323961e-9
1266 + tn[4] * 1.34942187e-12;
1267 smh[23] = 4.70720924 - ti * 12841.6265 + tn[0] * 4.30646568 - tn[1] *
1268 .00209329446 + tn[2] * 8.28571345e-6 - tn[3] * 4.99272172e-9
1269 + tn[4] * 1.15254502e-12;
1270 smh[24] = ti * 11522.2055 + 2.66682316 + tn[0] * 4.29142492 - tn[1] *
1271 .00275077135 + tn[2] * 9.99063813e-6 - tn[3] * 5.90388571e-9
1272 + tn[4] * 1.34342886e-12;
1273 smh[25] = 12.490417 - ti * 20059.449 + tn[0] * 2.2517214 + tn[1] *
1274 .0088275105 - tn[2] * 3.95485017e-6 + tn[3] * 1.43964658e-9 -
1275 tn[4] * 2.53324055e-13;
1276 smh[26] = ti * 7042.91804 + 12.215648 + tn[0] * 2.1358363 + tn[1] *
1277 .00905943605 - tn[2] * 2.89912457e-6 + tn[3] * 7.7866464e-10
1278 - tn[4] * 1.00728807e-13;
1279 smh[27] = 9.5714535 - ti * 1521.4766 + tn[0] * 3.4090624 + tn[1] *
1280 .005369287 + tn[2] * 3.1524875e-7 + tn[3] * 5.96548592e-10 +
1281 tn[4] * 1.43369255e-13;
1282 smh[28] = 16.14534 - ti * 1074.826 + tn[0] * 1.493307 + tn[1] *
1283 .01046259 + tn[2] * 7.47799e-7 - tn[3] * 1.39076e-9 + tn[4] *
1284 3.579073e-13;
1285 }
1286
1287 return 0;
1288} /* rdsmh_ */
1289
1290/* C */
1291/* ----------------------------------------------------------------------C */
1292/* C */
1293/* Subroutine */ int ratx_(double *t, double *alogt, double *c__,
1294 double *rf, double *rb, double *rklow)
1295{
1296 /* Initialized data */
1297
1298 static double small = 1e-50;
1299
1300 /* System generated locals */
1301 double d__1;
1302
1303 /* Builtin functions */
1304 //dol double log10(double *), exp(double), pow(double *,
1305 //dol double *);
1306
1307 /* Local variables */
1308 static int k;
1309 static double fc, pr, xn, ctb[167], flog, pcor, ctot, fclog, fcent,
1310 prlog, cprlog;
1311
1312
1313/* *****precision > double */
1314/* *****END precision > double */
1315
1316 /* Parameter adjustments */
1317 --rklow;
1318 --rb;
1319 --rf;
1320 --c__;
1321
1322 /* Function Body */
1323
1324/* third-body reactions */
1325
1326 ctot = 0.f;
1327 for (k = 1; k <= 19; ++k) {
1328 ctot += c__[k];
1329 }
1330
1331 ctb[0] = ctot + c__[1] * 1.4 + c__[6] * 14.4 + c__[10] + c__[11] * .75 +
1332 c__[12] * 2.6 + c__[16] * 2. + c__[15] * 2. + c__[18] * 3.;
1333 ctb[1] = ctot + c__[1] + c__[6] * 5. + c__[10] + c__[11] * .5 + c__[12] +
1334 c__[16] * 2. + c__[15] * 2. + c__[18] * 3.;
1335 ctb[11] = ctot + c__[1] + c__[4] * 5. + c__[6] * 5. + c__[10] + c__[11] *
1336 .5 + c__[12] * 2.5 + c__[16] * 2. + c__[15] * 2. + c__[18] * 3.;
1337 ctb[27] = ctot - c__[4] - c__[6] - c__[11] * .25 + c__[12] * .5 + c__[16]
1338 * .5 - c__[19] + c__[15] * 2. + c__[18] * 3.;
1339 ctb[32] = ctot - c__[1] - c__[6] + c__[10] - c__[12] + c__[16] * 2. + c__[
1340 15] * 2. + c__[18] * 3.;
1341 ctb[36] = ctot - c__[1] * .27 + c__[6] * 2.65 + c__[10] + c__[16] * 2. +
1342 c__[15] * 2. + c__[18] * 3.;
1343 ctb[43] = ctb[1];
1344 ctb[45] = ctot + c__[1] + c__[6] * 5. + c__[10] * 2. + c__[11] * .5 + c__[
1345 12] + c__[16] * 2. + c__[15] * 2. + c__[18] * 3.;
1346 ctb[47] = ctb[1];
1347 ctb[49] = ctb[1];
1348 ctb[54] = ctb[1];
1349 ctb[55] = ctb[1];
1350 ctb[57] = ctb[1];
1351 ctb[59] = ctb[1];
1352 ctb[65] = ctb[1];
1353 ctb[67] = ctb[1];
1354 ctb[105] = ctb[1];
1355 ctb[114] = ctb[1];
1356 ctb[130] = ctb[1];
1357 ctb[137] = ctot + c__[1] - c__[6] + c__[10] + c__[11] * .5 + c__[12] +
1358 c__[16] * 2.;
1359 ctb[141] = ctb[1];
1360 ctb[149] = ctb[1];
1361 ctb[156] = ctb[1];
1362 ctb[166] = ctb[1];
1363
1364/* If fall-off (pressure correction): */
1365
1366
1367 pr = rklow[1] * ctb[11] / rf[12];
1368 pcor = pr / (pr + 1.f);
1369 rf[12] *= pcor;
1370 rb[12] *= pcor;
1371
1372 pr = rklow[2] * ctb[43] / rf[44];
1373 pcor = pr / (pr + 1.f);
1374 d__1 = max(pr,small);
1375 prlog = log10(d__1);
1376 fcent = exp(-(*t) / 91.) * .438 + exp(-(*t) / 5836.) * .562 + exp(-8552. /
1377 *t);
1378 d__1 = max(fcent,small);
1379 fclog = log10(d__1);
1380 xn = .75f - fclog * 1.27f;
1381 cprlog = prlog - (fclog * .67f + .4f);
1382/* Computing 2nd power */
1383 d__1 = cprlog / (xn - cprlog * .14f);
1384 flog = fclog / (d__1 * d__1 + 1.f);
1385 fc = pow(c_b5, flog);
1386 pcor = fc * pcor;
1387 rf[44] *= pcor;
1388 rb[44] *= pcor;
1389
1390 pr = rklow[3] * ctb[45] / rf[46];
1391 pcor = pr / (pr + 1.f);
1392 d__1 = max(pr,small);
1393 prlog = log10(d__1);
1394 fcent = exp(-(*t) / 74.) * .217 + exp(-(*t) / 2941.) * .783 + exp(-6964. /
1395 *t);
1396 d__1 = max(fcent,small);
1397 fclog = log10(d__1);
1398 xn = .75f - fclog * 1.27f;
1399 cprlog = prlog - (fclog * .67f + .4f);
1400/* Computing 2nd power */
1401 d__1 = cprlog / (xn - cprlog * .14f);
1402 flog = fclog / (d__1 * d__1 + 1.f);
1403 fc = pow(c_b5, flog);
1404 pcor = fc * pcor;
1405 rf[46] *= pcor;
1406 rb[46] *= pcor;
1407
1408 pr = rklow[4] * ctb[47] / rf[48];
1409 pcor = pr / (pr + 1.f);
1410 d__1 = max(pr,small);
1411 prlog = log10(d__1);
1412 fcent = exp(-(*t) / 271.) * .2176 + exp(-(*t) / 2755.) * .7824 + exp(
1413 -6570. / *t);
1414 d__1 = max(fcent,small);
1415 fclog = log10(d__1);
1416 xn = .75f - fclog * 1.27f;
1417 cprlog = prlog - (fclog * .67f + .4f);
1418/* Computing 2nd power */
1419 d__1 = cprlog / (xn - cprlog * .14f);
1420 flog = fclog / (d__1 * d__1 + 1.f);
1421 fc = pow(c_b5, flog);
1422 pcor = fc * pcor;
1423 rf[48] *= pcor;
1424 rb[48] *= pcor;
1425
1426 pr = rklow[5] * ctb[49] / rf[50];
1427 pcor = pr / (pr + 1.f);
1428 d__1 = max(pr,small);
1429 prlog = log10(d__1);
1430 fcent = exp(-(*t) / 94.) * .242 + exp(-(*t) / 1555.) * .758 + exp(-4200. /
1431 *t);
1432 d__1 = max(fcent,small);
1433 fclog = log10(d__1);
1434 xn = .75f - fclog * 1.27f;
1435 cprlog = prlog - (fclog * .67f + .4f);
1436/* Computing 2nd power */
1437 d__1 = cprlog / (xn - cprlog * .14f);
1438 flog = fclog / (d__1 * d__1 + 1.f);
1439 fc = pow(c_b5, flog);
1440 pcor = fc * pcor;
1441 rf[50] *= pcor;
1442 rb[50] *= pcor;
1443
1444 pr = rklow[6] * ctb[54] / rf[55];
1445 pcor = pr / (pr + 1.f);
1446 d__1 = max(pr,small);
1447 prlog = log10(d__1);
1448 fcent = exp(-(*t) / 98.5) * .2493 + exp(-(*t) / 1302.) * .7507 + exp(
1449 -4167. / *t);
1450 d__1 = max(fcent,small);
1451 fclog = log10(d__1);
1452 xn = .75f - fclog * 1.27f;
1453 cprlog = prlog - (fclog * .67f + .4f);
1454/* Computing 2nd power */
1455 d__1 = cprlog / (xn - cprlog * .14f);
1456 flog = fclog / (d__1 * d__1 + 1.f);
1457 fc = pow(c_b5, flog);
1458 pcor = fc * pcor;
1459 rf[55] *= pcor;
1460 rb[55] *= pcor;
1461
1462 pr = rklow[7] * ctb[55] / rf[56];
1463 pcor = pr / (pr + 1.f);
1464 d__1 = max(pr,small);
1465 prlog = log10(d__1);
1466 fcent = exp(-(*t) / 207.5) * .218 + exp(-(*t) / 2663.) * .782 + exp(
1467 -6095. / *t);
1468 d__1 = max(fcent,small);
1469 fclog = log10(d__1);
1470 xn = .75f - fclog * 1.27f;
1471 cprlog = prlog - (fclog * .67f + .4f);
1472/* Computing 2nd power */
1473 d__1 = cprlog / (xn - cprlog * .14f);
1474 flog = fclog / (d__1 * d__1 + 1.f);
1475 fc = pow(c_b5, flog);
1476 pcor = fc * pcor;
1477 rf[56] *= pcor;
1478 rb[56] *= pcor;
1479
1480 pr = rklow[8] * ctb[57] / rf[58];
1481 pcor = pr / (pr + 1.f);
1482 d__1 = max(pr,small);
1483 prlog = log10(d__1);
1484 fcent = exp(-(*t) / 210.) * .0247 + exp(-(*t) / 984.) * .9753 + exp(
1485 -4374. / *t);
1486 d__1 = max(fcent,small);
1487 fclog = log10(d__1);
1488 xn = .75f - fclog * 1.27f;
1489 cprlog = prlog - (fclog * .67f + .4f);
1490/* Computing 2nd power */
1491 d__1 = cprlog / (xn - cprlog * .14f);
1492 flog = fclog / (d__1 * d__1 + 1.f);
1493 fc = pow(c_b5, flog);
1494 pcor = fc * pcor;
1495 rf[58] *= pcor;
1496 rb[58] *= pcor;
1497
1498 pr = rklow[9] * ctb[59] / rf[60];
1499 pcor = pr / (pr + 1.f);
1500 d__1 = max(pr,small);
1501 prlog = log10(d__1);
1502 fcent = exp(-(*t) / 125.) * .1578 + exp(-(*t) / 2219.) * .8422 + exp(
1503 -6882. / *t);
1504 d__1 = max(fcent,small);
1505 fclog = log10(d__1);
1506 xn = .75f - fclog * 1.27f;
1507 cprlog = prlog - (fclog * .67f + .4f);
1508/* Computing 2nd power */
1509 d__1 = cprlog / (xn - cprlog * .14f);
1510 flog = fclog / (d__1 * d__1 + 1.f);
1511 fc = pow(c_b5, flog);
1512 pcor = fc * pcor;
1513 rf[60] *= pcor;
1514 rb[60] *= pcor;
1515
1516 pr = rklow[10] * ctb[65] / rf[66];
1517 pcor = pr / (pr + 1.f);
1518 d__1 = max(pr,small);
1519 prlog = log10(d__1);
1520 fcent = exp(-(*t) / 197.) * .068 + exp(-(*t) / 1540.) * .932 + exp(
1521 -10300. / *t);
1522 d__1 = max(fcent,small);
1523 fclog = log10(d__1);
1524 xn = .75f - fclog * 1.27f;
1525 cprlog = prlog - (fclog * .67f + .4f);
1526/* Computing 2nd power */
1527 d__1 = cprlog / (xn - cprlog * .14f);
1528 flog = fclog / (d__1 * d__1 + 1.f);
1529 fc = pow(c_b5, flog);
1530 pcor = fc * pcor;
1531 rf[66] *= pcor;
1532 rb[66] *= pcor;
1533
1534 pr = rklow[11] * ctb[67] / rf[68];
1535 pcor = pr / (pr + 1.f);
1536 d__1 = max(pr,small);
1537 prlog = log10(d__1);
1538 fcent = exp(-(*t) / 94.) * .2654 + exp(-(*t) / 1756.) * .7346 + exp(
1539 -5182. / *t);
1540 d__1 = max(fcent,small);
1541 fclog = log10(d__1);
1542 xn = .75f - fclog * 1.27f;
1543 cprlog = prlog - (fclog * .67f + .4f);
1544/* Computing 2nd power */
1545 d__1 = cprlog / (xn - cprlog * .14f);
1546 flog = fclog / (d__1 * d__1 + 1.f);
1547 fc = pow(c_b5, flog);
1548 pcor = fc * pcor;
1549 rf[68] *= pcor;
1550 rb[68] *= pcor;
1551
1552 pr = rklow[12] * ctb[105] / rf[106];
1553 pcor = pr / (pr + 1.f);
1554 d__1 = max(pr,small);
1555 prlog = log10(d__1);
1556 fcent = exp(-(*t) / 237.) * .4243 + exp(-(*t) / 1652.) * .5757 + exp(
1557 -5069. / *t);
1558 d__1 = max(fcent,small);
1559 fclog = log10(d__1);
1560 xn = .75f - fclog * 1.27f;
1561 cprlog = prlog - (fclog * .67f + .4f);
1562/* Computing 2nd power */
1563 d__1 = cprlog / (xn - cprlog * .14f);
1564 flog = fclog / (d__1 * d__1 + 1.f);
1565 fc = pow(c_b5, flog);
1566 pcor = fc * pcor;
1567 rf[106] *= pcor;
1568 rb[106] *= pcor;
1569
1570 pr = rklow[13] * ctb[114] / rf[115];
1571 pcor = pr / (pr + 1.f);
1572 d__1 = max(pr,small);
1573 prlog = log10(d__1);
1574 fcent = exp(-(*t) / 275.) * .4093 + exp(-(*t) / 1226.) * .5907 + exp(
1575 -5185. / *t);
1576 d__1 = max(fcent,small);
1577 fclog = log10(d__1);
1578 xn = .75f - fclog * 1.27f;
1579 cprlog = prlog - (fclog * .67f + .4f);
1580/* Computing 2nd power */
1581 d__1 = cprlog / (xn - cprlog * .14f);
1582 flog = fclog / (d__1 * d__1 + 1.f);
1583 fc = pow(c_b5, flog);
1584 pcor = fc * pcor;
1585 rf[115] *= pcor;
1586 rb[115] *= pcor;
1587
1588 pr = rklow[14] * ctb[130] / rf[131];
1589 pcor = pr / (pr + 1.f);
1590 d__1 = max(pr,small);
1591 prlog = log10(d__1);
1592 fcent = exp(-(*t) / 73.2) * .381 + exp(-(*t) / 1180.) * .619 + exp(-9999.
1593 / *t);
1594 d__1 = max(fcent,small);
1595 fclog = log10(d__1);
1596 xn = .75f - fclog * 1.27f;
1597 cprlog = prlog - (fclog * .67f + .4f);
1598/* Computing 2nd power */
1599 d__1 = cprlog / (xn - cprlog * .14f);
1600 flog = fclog / (d__1 * d__1 + 1.f);
1601 fc = pow(c_b5, flog);
1602 pcor = fc * pcor;
1603 rf[131] *= pcor;
1604 rb[131] *= pcor;
1605
1606 pr = rklow[15] * ctb[141] / rf[142];
1607 pcor = pr / (pr + 1.f);
1608 d__1 = max(pr,small);
1609 prlog = log10(d__1);
1610 fcent = exp(-(*t) / 180.) * .2655 + exp(-(*t) / 1035.) * .7345 + exp(
1611 -5417. / *t);
1612 d__1 = max(fcent,small);
1613 fclog = log10(d__1);
1614 xn = .75f - fclog * 1.27f;
1615 cprlog = prlog - (fclog * .67f + .4f);
1616/* Computing 2nd power */
1617 d__1 = cprlog / (xn - cprlog * .14f);
1618 flog = fclog / (d__1 * d__1 + 1.f);
1619 fc = pow(c_b5, flog);
1620 pcor = fc * pcor;
1621 rf[142] *= pcor;
1622 rb[142] *= pcor;
1623
1624 pr = rklow[16] * ctb[149] / rf[150];
1625 pcor = pr / (pr + 1.f);
1626 d__1 = max(pr,small);
1627 prlog = log10(d__1);
1628 fcent = exp(-(*t) / 122.) * .422 + exp(-(*t) / 2535.) * .578 + exp(-9365.
1629 / *t);
1630 d__1 = max(fcent,small);
1631 fclog = log10(d__1);
1632 xn = .75f - fclog * 1.27f;
1633 cprlog = prlog - (fclog * .67f + .4f);
1634/* Computing 2nd power */
1635 d__1 = cprlog / (xn - cprlog * .14f);
1636 flog = fclog / (d__1 * d__1 + 1.f);
1637 fc = pow(c_b5, flog);
1638 pcor = fc * pcor;
1639 rf[150] *= pcor;
1640 rb[150] *= pcor;
1641
1642 pr = rklow[17] * ctb[156] / rf[157];
1643 pcor = pr / (pr + 1.f);
1644 d__1 = max(pr,small);
1645 prlog = log10(d__1);
1646 fcent = exp(-(*t) / 201.) * .535 + exp(-(*t) / 1773.) * .465 + exp(-5333.
1647 / *t);
1648 d__1 = max(fcent,small);
1649 fclog = log10(d__1);
1650 xn = .75f - fclog * 1.27f;
1651 cprlog = prlog - (fclog * .67f + .4f);
1652/* Computing 2nd power */
1653 d__1 = cprlog / (xn - cprlog * .14f);
1654 flog = fclog / (d__1 * d__1 + 1.f);
1655 fc = pow(c_b5, flog);
1656 pcor = fc * pcor;
1657 rf[157] *= pcor;
1658 rb[157] *= pcor;
1659
1660 pr = rklow[18] * ctb[166] / rf[167];
1661 pcor = pr / (pr + 1.f);
1662 d__1 = max(pr,small);
1663 prlog = log10(d__1);
1664 fcent = exp(-(*t) / 1340.6) * .825 + exp(-(*t) / 6e4) * .175 + exp(
1665 -10139.8 / *t);
1666 d__1 = max(fcent,small);
1667 fclog = log10(d__1);
1668 xn = .75f - fclog * 1.27f;
1669 cprlog = prlog - (fclog * .67f + .4f);
1670/* Computing 2nd power */
1671 d__1 = cprlog / (xn - cprlog * .14f);
1672 flog = fclog / (d__1 * d__1 + 1.f);
1673 fc = pow(c_b5, flog);
1674 pcor = fc * pcor;
1675 rf[167] *= pcor;
1676 rb[167] *= pcor;
1677
1678 rf[1] = rf[1] * ctb[0] * c__[3] * c__[3];
1679 rf[2] = rf[2] * ctb[1] * c__[3] * c__[2];
1680 rf[3] = rf[3] * c__[3] * c__[1];
1681 rf[4] = rf[4] * c__[3] * c__[7];
1682 rf[5] = rf[5] * c__[3] * c__[8];
1683 rf[6] *= c__[3];
1684 rf[7] *= c__[3];
1685 rf[8] *= c__[3];
1686 rf[9] *= c__[3];
1687 rf[10] = rf[10] * c__[3] * c__[9];
1688 rf[11] = rf[11] * c__[3] * c__[10];
1689 rf[12] = rf[12] * c__[3] * c__[11];
1690 rf[13] *= c__[3];
1691 rf[14] *= c__[3];
1692 rf[15] = rf[15] * c__[3] * c__[13];
1693 rf[16] *= c__[3];
1694 rf[17] = rf[17] * c__[3] * c__[14];
1695 rf[18] = rf[18] * c__[3] * c__[14];
1696 rf[19] *= c__[3];
1697 rf[20] = rf[20] * c__[3] * c__[15];
1698 rf[21] *= c__[3];
1699 rf[22] = rf[22] * c__[3] * c__[16];
1700 rf[23] *= c__[3];
1701 rf[24] = rf[24] * c__[3] * c__[17];
1702 rf[25] = rf[25] * c__[3] * c__[17];
1703 rf[26] = rf[26] * c__[4] * c__[11];
1704 rf[27] = rf[27] * c__[4] * c__[13];
1705 rf[28] = rf[28] * ctb[27] * c__[2] * c__[4];
1706 rf[29] = rf[29] * c__[2] * c__[4] * c__[4];
1707 rf[30] = rf[30] * c__[2] * c__[4] * c__[6];
1708 rf[31] = rf[31] * c__[2] * c__[4] * c__[19];
1709 rf[32] = rf[32] * c__[2] * c__[4];
1710 rf[33] = rf[33] * ctb[32] * c__[2] * c__[2];
1711 rf[34] = rf[34] * c__[2] * c__[2] * c__[1];
1712 rf[35] = rf[35] * c__[2] * c__[2] * c__[6];
1713 rf[36] = rf[36] * c__[2] * c__[2] * c__[12];
1714 rf[37] = rf[37] * ctb[36] * c__[2] * c__[5];
1715 rf[38] = rf[38] * c__[2] * c__[7];
1716 rf[39] = rf[39] * c__[2] * c__[7];
1717 rf[40] = rf[40] * c__[2] * c__[7];
1718 rf[41] = rf[41] * c__[2] * c__[8];
1719 rf[42] = rf[42] * c__[2] * c__[8];
1720 rf[43] *= c__[2];
1721 rf[44] *= c__[2];
1722 rf[45] *= c__[2];
1723 rf[46] = rf[46] * c__[2] * c__[9];
1724 rf[47] = rf[47] * c__[2] * c__[10];
1725 rf[48] *= c__[2];
1726 rf[49] *= c__[2];
1727 rf[50] = rf[50] * c__[2] * c__[13];
1728 rf[51] = rf[51] * c__[2] * c__[13];
1729 rf[52] *= c__[2];
1730 rf[53] *= c__[2];
1731 rf[54] *= c__[2];
1732 rf[55] = rf[55] * c__[2] * c__[14];
1733 rf[56] *= c__[2];
1734 rf[57] *= c__[2];
1735 rf[58] = rf[58] * c__[2] * c__[15];
1736 rf[59] = rf[59] * c__[2] * c__[15];
1737 rf[60] *= c__[2];
1738 rf[61] *= c__[2];
1739 rf[62] = rf[62] * c__[2] * c__[16];
1740 rf[63] *= c__[2];
1741 rf[64] = rf[64] * c__[2] * c__[17];
1742 rf[65] = rf[65] * c__[2] * c__[17];
1743 rf[66] = rf[66] * c__[1] * c__[11];
1744 rf[67] = rf[67] * c__[5] * c__[1];
1745 rf[68] = rf[68] * c__[5] * c__[5];
1746 rf[69] = rf[69] * c__[5] * c__[5];
1747 rf[70] = rf[70] * c__[5] * c__[7];
1748 rf[71] = rf[71] * c__[5] * c__[8];
1749 rf[72] = rf[72] * c__[5] * c__[8];
1750 rf[73] *= c__[5];
1751 rf[74] *= c__[5];
1752 rf[75] *= c__[5];
1753 rf[76] *= c__[5];
1754 rf[77] *= c__[5];
1755 rf[78] = rf[78] * c__[5] * c__[9];
1756 rf[79] = rf[79] * c__[5] * c__[9];
1757 rf[80] = rf[80] * c__[5] * c__[10];
1758 rf[81] = rf[81] * c__[5] * c__[11];
1759 rf[82] *= c__[5];
1760 rf[83] = rf[83] * c__[5] * c__[13];
1761 rf[84] *= c__[5];
1762 rf[85] = rf[85] * c__[5] * c__[14];
1763 rf[86] = rf[86] * c__[5] * c__[14];
1764 rf[87] *= c__[5];
1765 rf[88] = rf[88] * c__[5] * c__[15];
1766 rf[89] = rf[89] * c__[5] * c__[16];
1767 rf[90] = rf[90] * c__[5] * c__[17];
1768 rf[91] = rf[91] * c__[7] * c__[7];
1769 rf[92] = rf[92] * c__[7] * c__[7];
1770 rf[93] *= c__[7];
1771 rf[94] = rf[94] * c__[7] * c__[9];
1772 rf[95] = rf[95] * c__[7] * c__[9];
1773 rf[96] = rf[96] * c__[7] * c__[11];
1774 rf[97] = rf[97] * c__[7] * c__[13];
1775 rf[98] *= c__[4];
1776 rf[99] *= c__[9];
1777 rf[100] *= c__[4];
1778 rf[101] *= c__[1];
1779 rf[102] *= c__[6];
1780 rf[105] *= c__[10];
1781 rf[106] *= c__[11];
1782 rf[107] *= c__[12];
1783 rf[108] *= c__[13];
1784 rf[110] *= c__[4];
1785 rf[111] *= c__[1];
1786 rf[113] *= c__[9];
1787 rf[114] *= c__[10];
1788 rf[115] *= c__[11];
1789 rf[117] *= c__[19];
1790 rf[118] *= c__[4];
1791 rf[119] *= c__[4];
1792 rf[120] *= c__[1];
1793 rf[121] *= c__[6];
1794 rf[122] *= c__[9];
1795 rf[123] *= c__[10];
1796 rf[124] *= c__[11];
1797 rf[125] *= c__[12];
1798 rf[126] *= c__[12];
1799 rf[128] = rf[128] * c__[9] * c__[4];
1800 rf[129] = rf[129] * c__[9] * c__[4];
1801 rf[130] = rf[130] * c__[9] * c__[8];
1802 rf[131] = rf[131] * c__[9] * c__[9];
1803 rf[132] = rf[132] * c__[9] * c__[9];
1804 rf[133] *= c__[9];
1805 rf[134] = rf[134] * c__[9] * c__[13];
1806 rf[135] = rf[135] * c__[9] * c__[15];
1807 rf[136] = rf[136] * c__[9] * c__[16];
1808 rf[137] *= c__[6];
1809 rf[138] *= ctb[137];
1810 rf[139] *= c__[4];
1811 rf[140] *= c__[4];
1812 rf[141] *= c__[4];
1813 rf[142] *= c__[15];
1814 rf[143] *= c__[4];
1815 rf[144] *= c__[4];
1816 rf[146] = rf[146] * c__[3] * c__[9];
1817 rf[147] = rf[147] * c__[3] * c__[15];
1818 rf[148] = rf[148] * c__[5] * c__[7];
1819 rf[149] = rf[149] * c__[5] * c__[9];
1820 rf[150] *= c__[1];
1821 rf[151] *= c__[4];
1822 rf[152] *= c__[4];
1823 rf[154] *= c__[6];
1824 rf[155] *= c__[4];
1825 rf[156] *= c__[4];
1826 rf[157] = rf[157] * c__[2] * c__[17];
1827 rf[159] *= c__[4];
1828 rf[162] *= c__[2];
1829 rf[163] *= c__[5];
1830 rf[164] = rf[164] * c__[18] * c__[2];
1831 rf[165] = rf[165] * c__[18] * c__[3];
1832 rf[166] = rf[166] * c__[18] * c__[3];
1833 rf[167] *= c__[9];
1834 rb[1] = rb[1] * ctb[0] * c__[4];
1835 rb[2] = rb[2] * ctb[1] * c__[5];
1836 rb[3] = rb[3] * c__[2] * c__[5];
1837 rb[4] = rb[4] * c__[5] * c__[4];
1838 rb[5] = rb[5] * c__[5] * c__[7];
1839 rb[6] = rb[6] * c__[2] * c__[11];
1840 rb[7] *= c__[2];
1841 rb[8] = rb[8] * c__[1] * c__[11];
1842 rb[9] *= c__[2];
1843 rb[10] = rb[10] * c__[2] * c__[13];
1844 rb[11] = rb[11] * c__[5] * c__[9];
1845 rb[12] *= c__[12];
1846 rb[13] = rb[13] * c__[5] * c__[11];
1847 rb[14] = rb[14] * c__[2] * c__[12];
1848 rb[15] *= c__[5];
1849 rb[16] = rb[16] * c__[5] * c__[13];
1850 rb[17] *= c__[2];
1851 rb[18] *= c__[11];
1852 rb[19] = rb[19] * c__[2] * c__[17];
1853 rb[20] *= c__[9];
1854 rb[21] = rb[21] * c__[9] * c__[13];
1855 rb[22] *= c__[5];
1856 rb[23] = rb[23] * c__[2] * c__[11] * c__[11];
1857 rb[24] *= c__[5];
1858 rb[25] *= c__[12];
1859 rb[26] = rb[26] * c__[3] * c__[12];
1860 rb[27] *= c__[7];
1861 rb[28] = rb[28] * ctb[27] * c__[7];
1862 rb[29] = rb[29] * c__[7] * c__[4];
1863 rb[30] = rb[30] * c__[7] * c__[6];
1864 rb[31] = rb[31] * c__[7] * c__[19];
1865 rb[32] = rb[32] * c__[3] * c__[5];
1866 rb[33] = rb[33] * ctb[32] * c__[1];
1867 rb[34] = rb[34] * c__[1] * c__[1];
1868 rb[35] = rb[35] * c__[1] * c__[6];
1869 rb[36] = rb[36] * c__[1] * c__[12];
1870 rb[37] = rb[37] * ctb[36] * c__[6];
1871 rb[38] = rb[38] * c__[3] * c__[6];
1872 rb[39] = rb[39] * c__[4] * c__[1];
1873 rb[40] = rb[40] * c__[5] * c__[5];
1874 rb[41] = rb[41] * c__[7] * c__[1];
1875 rb[42] = rb[42] * c__[5] * c__[6];
1876 rb[43] *= c__[1];
1877 rb[44] *= c__[9];
1878 rb[45] *= c__[1];
1879 rb[46] *= c__[10];
1880 rb[47] = rb[47] * c__[9] * c__[1];
1881 rb[48] *= c__[13];
1882 rb[49] = rb[49] * c__[1] * c__[11];
1883 rb[51] *= c__[1];
1884 rb[52] = rb[52] * c__[1] * c__[13];
1885 rb[53] = rb[53] * c__[5] * c__[9];
1886 rb[54] *= c__[6];
1887 rb[56] *= c__[15];
1888 rb[57] = rb[57] * c__[1] * c__[14];
1889 rb[59] *= c__[1];
1890 rb[60] *= c__[16];
1891 rb[61] = rb[61] * c__[1] * c__[15];
1892 rb[62] *= c__[1];
1893 rb[63] *= c__[11];
1894 rb[64] *= c__[1];
1895 rb[65] = rb[65] * c__[9] * c__[11];
1896 rb[66] *= c__[13];
1897 rb[67] = rb[67] * c__[2] * c__[6];
1898 rb[68] *= c__[8];
1899 rb[69] = rb[69] * c__[3] * c__[6];
1900 rb[70] = rb[70] * c__[4] * c__[6];
1901 rb[71] = rb[71] * c__[7] * c__[6];
1902 rb[72] = rb[72] * c__[7] * c__[6];
1903 rb[73] = rb[73] * c__[2] * c__[11];
1904 rb[74] *= c__[2];
1905 rb[75] = rb[75] * c__[2] * c__[13];
1906 rb[76] *= c__[6];
1907 rb[77] = rb[77] * c__[2] * c__[13];
1908 rb[78] *= c__[6];
1909 rb[79] *= c__[6];
1910 rb[80] = rb[80] * c__[9] * c__[6];
1911 rb[81] = rb[81] * c__[2] * c__[12];
1912 rb[82] = rb[82] * c__[6] * c__[11];
1913 rb[83] *= c__[6];
1914 rb[84] = rb[84] * c__[6] * c__[13];
1915 rb[85] = rb[85] * c__[2] * c__[17];
1916 rb[86] = rb[86] * c__[9] * c__[11];
1917 rb[87] = rb[87] * c__[6] * c__[14];
1918 rb[88] *= c__[6];
1919 rb[89] *= c__[6];
1920 rb[90] *= c__[6];
1921 rb[91] = rb[91] * c__[4] * c__[8];
1922 rb[92] = rb[92] * c__[4] * c__[8];
1923 rb[93] = rb[93] * c__[5] * c__[13];
1924 rb[94] = rb[94] * c__[4] * c__[10];
1925 rb[95] *= c__[5];
1926 rb[96] = rb[96] * c__[5] * c__[12];
1927 rb[97] *= c__[8];
1928 rb[98] = rb[98] * c__[3] * c__[11];
1929 rb[99] = rb[99] * c__[2] * c__[14];
1930 rb[100] *= c__[3];
1931 rb[101] *= c__[2];
1932 rb[102] = rb[102] * c__[2] * c__[13];
1933 rb[103] = rb[103] * c__[2] * c__[14];
1934 rb[105] = rb[105] * c__[2] * c__[15];
1935 rb[107] *= c__[11];
1936 rb[108] = rb[108] * c__[2] * c__[17];
1937 rb[109] = rb[109] * c__[11] * c__[14];
1938 rb[111] = rb[111] * c__[2] * c__[9];
1939 rb[112] = rb[112] * c__[1] * c__[14];
1940 rb[113] = rb[113] * c__[2] * c__[15];
1941 rb[114] = rb[114] * c__[9] * c__[9];
1942 rb[115] *= c__[17];
1943 rb[117] *= c__[19];
1944 rb[118] = rb[118] * c__[2] * c__[5] * c__[11];
1945 rb[119] = rb[119] * c__[11] * c__[6];
1946 rb[120] = rb[120] * c__[9] * c__[2];
1947 rb[121] *= c__[6];
1948 rb[122] = rb[122] * c__[2] * c__[15];
1949 rb[123] = rb[123] * c__[9] * c__[9];
1950 rb[124] *= c__[11];
1951 rb[125] *= c__[12];
1952 rb[126] = rb[126] * c__[11] * c__[13];
1953 rb[128] *= c__[3];
1954 rb[129] = rb[129] * c__[5] * c__[13];
1955 rb[130] = rb[130] * c__[7] * c__[10];
1956 rb[131] *= c__[16];
1957 rb[132] *= c__[2];
1958 rb[133] = rb[133] * c__[10] * c__[11];
1959 rb[134] *= c__[10];
1960 rb[135] *= c__[10];
1961 rb[136] *= c__[10];
1962 rb[137] = rb[137] * c__[2] * c__[11] * c__[6];
1963 rb[138] = rb[138] * ctb[137] * c__[2] * c__[11];
1964 rb[139] = rb[139] * c__[7] * c__[11];
1965 rb[140] = rb[140] * c__[7] * c__[13];
1966 rb[142] = rb[142] * c__[1] * c__[14];
1967 rb[143] = rb[143] * c__[7] * c__[15];
1968 rb[144] = rb[144] * c__[5] * c__[11] * c__[11];
1969 rb[145] = rb[145] * c__[11] * c__[11] * c__[14];
1970 rb[147] *= c__[2];
1971 rb[148] = rb[148] * c__[4] * c__[6];
1972 rb[150] *= c__[9];
1973 rb[152] = rb[152] * c__[3] * c__[13];
1974 rb[156] = rb[156] * c__[7] * c__[14];
1975 rb[161] *= c__[9];
1976 rb[162] = rb[162] * c__[17] * c__[1];
1977 rb[163] = rb[163] * c__[6] * c__[17];
1978 rb[164] = rb[164] * c__[15] * c__[9];
1979 rb[165] = rb[165] * c__[17] * c__[9] * c__[2];
1980 rb[167] *= c__[18];
1981
1982 return 0;
1983} /* ratx_ */
1984
int rdsmh_(double *t, double *tlog, double *smh)
Definition c2h4red.cc:1067
int ratx_(double *t, double *alogt, double *c__, double *rf, double *rb, double *rklow)
Definition c2h4red.cc:1293
int ratt_(double *t, double *alogt, double *rf, double *rb, double *rklow)
Definition c2h4red.cc:514
void getProblemSpecificRR(double rho, double temp, double pres, double *yi, double *rrsp)
Definition c2h4red.cc:47