1 |
REQUIRE "atoms.a4l"; |
2 |
(* => atoms.a4l, measures.a4l, system.a4l, basemodel.a4l *) |
3 |
PROVIDE "kinetics.a4l"; |
4 |
|
5 |
(* |
6 |
* kinetics.a4l |
7 |
* by Duncan Coffey |
8 |
* Part of the ASCEND Library |
9 |
* $Date: 1998/06/17 19:10:29 $ |
10 |
* $Revision: 1.3 $ |
11 |
* $Author: mthomas $ |
12 |
* $Source: /afs/cs.cmu.edu/project/ascend/Repository/models/kinetics.a4l,v $ |
13 |
* |
14 |
* This file is part of the ASCEND Modeling Library. |
15 |
* |
16 |
* Copyright (C) 1998 Duncan Coffey |
17 |
* |
18 |
* The ASCEND Modeling Library is free software; you can redistribute |
19 |
* it and/or modify it under the terms of the GNU General Public |
20 |
* License as published by the Free Software Foundation; either |
21 |
* version 2 of the License, or (at your option) any later version. |
22 |
* |
23 |
* The ASCEND Modeling Library is distributed in hope that it |
24 |
* will be useful, but WITHOUT ANY WARRANTY; without even the implied |
25 |
* warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. |
26 |
* See the GNU General Public License for more details. |
27 |
* |
28 |
* You should have received a copy of the GNU General Public License |
29 |
* along with the program; if not, write to the Free Software |
30 |
* Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139 USA. Check |
31 |
* the file named COPYING. |
32 |
*) |
33 |
|
34 |
|
35 |
(* ************************************************************************ *) |
36 |
(* ******************** Kinetics ******************** *) |
37 |
(* ************************************************************************ *) |
38 |
MODEL base_kinetics( |
39 |
components IS_A set OF symbol_constant; |
40 |
nr IS_A set OF symbol_constant; |
41 |
Kr[nr] IS_A constant; |
42 |
active_energy[nr] IS_A constant; |
43 |
reac_T WILL_BE temperature; |
44 |
species[nr] IS_A set OF symbol_constant; |
45 |
nu[components][nr] WILL_BE integer_constant; |
46 |
conc[components] WILL_BE molar_density; |
47 |
)WHERE( |
48 |
|
49 |
)REFINES cmumodel; |
50 |
NOTES |
51 |
'purpose' SELF { |
52 |
This MODEL computes the production rates per volume of the reactant |
53 |
speices given. There are three possible models to use. |
54 |
|
55 |
base_kinetics - No rate law is constructed and is expected to be given |
56 |
by the user. Production rates are then calculated from the given rate laws. |
57 |
|
58 |
nr - names of reactions |
59 |
Kr[nr] - rate constants |
60 |
species[nr] - set of components involved in each reaction |
61 |
nu[components][nr] - stoichiometry for each components and reaction |
62 |
conc[components] - Concentrations used in rate laws. |
63 |
|
64 |
element_kinetics - Here the rate laws are constructed from the |
65 |
stoichiometry (nu) given. |
66 |
|
67 |
specify_kinetics - Here the rate laws are constructed from the orders (order) given. |
68 |
|
69 |
} |
70 |
'developer-Duncan' SELF {From the element_kinetics and |
71 |
specify_kinetics most rate laws can be constructed. However in many cases |
72 |
there may need to be intermediate steps and it might just be easier to |
73 |
use base_kinetics. |
74 |
} |
75 |
END NOTES; |
76 |
|
77 |
rate[nr] IS_A conc_rate; |
78 |
reactions[components] IS_A set OF symbol_constant; |
79 |
production[components] IS_A conc_rate; |
80 |
Ftot_scale IS_A scaling_constant; |
81 |
R IS_A gas_constant; |
82 |
|
83 |
(* Define the production rate for a component *) |
84 |
FOR i IN components CREATE |
85 |
production[i] / Ftot_scale = SUM[nu[i][j]*rate[j] | j IN reactions[i]] / Ftot_scale; |
86 |
END FOR; |
87 |
|
88 |
(* define rxns where each component is present *) |
89 |
FOR i IN components CREATE |
90 |
reactions[i] :== [j IN nr | i IN species[j]]; |
91 |
END FOR; |
92 |
|
93 |
METHODS |
94 |
|
95 |
METHOD default_self; |
96 |
rate[nr].lower_bound :=-1e100 {mole/m^3/s}; |
97 |
production[components].lower_bound :=-1e100 {mole/m^3/s}; |
98 |
Ftot_scale :=1000; |
99 |
END default_self; |
100 |
|
101 |
METHOD default_all; |
102 |
RUN reac.default_self; |
103 |
END default_all; |
104 |
|
105 |
METHOD check_self; |
106 |
END check_self; |
107 |
|
108 |
METHOD check_all; |
109 |
RUN check_self; |
110 |
END check_all; |
111 |
|
112 |
METHOD bound_self; |
113 |
END bound_self; |
114 |
|
115 |
METHOD bound_all; |
116 |
RUN bound_self; |
117 |
END bound_all; |
118 |
|
119 |
METHOD scale_self; |
120 |
END scale_self; |
121 |
|
122 |
METHOD scale_all; |
123 |
RUN scale_self; |
124 |
END scale_all; |
125 |
|
126 |
METHOD seqmod; |
127 |
conc[components].fixed :=TRUE; |
128 |
reac_T.fixed :=TRUE; |
129 |
END seqmod; |
130 |
|
131 |
METHOD specify; |
132 |
RUN seqmod; |
133 |
END specify; |
134 |
END base_kinetics; |
135 |
|
136 |
|
137 |
MODEL element_kinetics( |
138 |
components IS_A set OF symbol_constant; |
139 |
nr IS_A set OF symbol_constant; |
140 |
Kr[nr] IS_A constant; |
141 |
active_energy[nr] IS_A constant; |
142 |
reac_T WILL_BE temperature; |
143 |
species[nr] IS_A set OF symbol_constant; |
144 |
nu[components][nr] WILL_BE integer_constant; |
145 |
conc[components] WILL_BE molar_density; |
146 |
)WHERE( |
147 |
|
148 |
)REFINES base_kinetics; |
149 |
|
150 |
(* define the rate equations for each component *) |
151 |
FOR j IN nr CREATE |
152 |
rate[j] / Ftot_scale = (Kr[j]*exp(-active_energy[j]/R/reac_T)*PROD[ |
153 |
PROD[ conc[i] | m IN [1..-(nu[i][j])]] | i IN species[j]]) / Ftot_scale; |
154 |
END FOR; |
155 |
|
156 |
|
157 |
|
158 |
METHODS |
159 |
|
160 |
METHOD default_self; |
161 |
rate[nr].lower_bound :=-1e100 {mole/m^3/s}; |
162 |
production[components].lower_bound :=-1e100 {mole/m^3/s}; |
163 |
Ftot_scale :=1000; |
164 |
END default_self; |
165 |
|
166 |
METHOD default_all; |
167 |
RUN reac.default_self; |
168 |
END default_all; |
169 |
|
170 |
METHOD check_self; |
171 |
END check_self; |
172 |
|
173 |
METHOD check_all; |
174 |
RUN check_self; |
175 |
END check_all; |
176 |
|
177 |
METHOD bound_self; |
178 |
END bound_self; |
179 |
|
180 |
METHOD bound_all; |
181 |
RUN bound_self; |
182 |
END bound_all; |
183 |
|
184 |
METHOD scale_self; |
185 |
END scale_self; |
186 |
|
187 |
METHOD scale_all; |
188 |
RUN scale_self; |
189 |
END scale_all; |
190 |
|
191 |
METHOD seqmod; |
192 |
conc[components].fixed :=TRUE; |
193 |
reac_T.fixed :=TRUE; |
194 |
END seqmod; |
195 |
|
196 |
METHOD specify; |
197 |
RUN seqmod; |
198 |
END specify; |
199 |
|
200 |
END element_kinetics; |
201 |
|
202 |
|
203 |
MODEL specify_kinetics( |
204 |
components IS_A set OF symbol_constant; |
205 |
nr IS_A set OF symbol_constant; |
206 |
Kr[nr] IS_A constant; |
207 |
active_energy[nr] IS_A constant; |
208 |
reac_T WILL_BE temperature; |
209 |
species[nr] IS_A set OF symbol_constant; |
210 |
nu[components][nr] WILL_BE integer_constant; |
211 |
conc[components] WILL_BE molar_density; |
212 |
order[components][nr] WILL_BE constant; |
213 |
)WHERE( |
214 |
|
215 |
)REFINES base_kinetics; |
216 |
|
217 |
|
218 |
(* define the rate equations for each component *) |
219 |
FOR j IN nr CREATE |
220 |
rate[j] / Ftot_scale = (Kr[j]*exp(-active_energy[j]/R/reac_T)* |
221 |
PROD[conc[i]^order[i][j] | i IN species[j]]) / Ftot_scale; |
222 |
END FOR; |
223 |
|
224 |
METHODS |
225 |
|
226 |
METHOD default_self; |
227 |
rate[nr].lower_bound :=-1e100 {mole/m^3/s}; |
228 |
production[components].lower_bound :=-1e100 {mole/m^3/s}; |
229 |
Ftot_scale :=1000; |
230 |
END default_self; |
231 |
|
232 |
METHOD seqmod; |
233 |
conc[components].fixed :=TRUE; |
234 |
reac_T.fixed :=TRUE; |
235 |
END seqmod; |
236 |
|
237 |
METHOD specify; |
238 |
RUN seqmod; |
239 |
END specify; |
240 |
|
241 |
END specify_kinetics; |
242 |
|
243 |
MODEL test_elm_kinetics; |
244 |
components IS_A set OF symbol_constant; |
245 |
nr IS_A set OF symbol_constant; |
246 |
Kr[nr] IS_A constant; |
247 |
active_energy[nr] IS_A constant; |
248 |
reac_T IS_A temperature; |
249 |
species[nr] IS_A set OF symbol_constant; |
250 |
nu[components][nr] IS_A integer_constant; |
251 |
conc[components] IS_A molar_density; |
252 |
|
253 |
|
254 |
components :==['methane','ethane','n_butane']; |
255 |
nr :==['reaction_1','reaction_2']; |
256 |
species['reaction_1'] :==['methane','ethane']; |
257 |
species['reaction_2'] :==['methane','ethane']; |
258 |
Kr['reaction_1'] :==2 {m^3/mol/s}; |
259 |
Kr['reaction_2'] :==3 {1/s}; |
260 |
active_energy['reaction_1'] :==70000 {J}; |
261 |
active_energy['reaction_2'] :==70000 {J}; |
262 |
|
263 |
nu['methane']['reaction_1'] :==-2; |
264 |
nu['ethane']['reaction_1'] :==1; |
265 |
nu['methane']['reaction_2'] :==2; |
266 |
nu['ethane']['reaction_2'] :==-1; |
267 |
(* define nu=0 FOR components not IN reactions*) |
268 |
FOR i IN nr CREATE |
269 |
FOR j IN components - [species[i]] CREATE |
270 |
nu[j][i] :==0; |
271 |
END FOR; |
272 |
END FOR; |
273 |
|
274 |
|
275 |
|
276 |
reac IS_A element_kinetics(components, nr, Kr, active_energy, |
277 |
reac_T, species, nu, conc); |
278 |
|
279 |
METHODS |
280 |
|
281 |
METHOD default_self; |
282 |
RUN reac.default_self; |
283 |
END default_self; |
284 |
|
285 |
METHOD default_all; |
286 |
RUN reac.default_self; |
287 |
END default_all; |
288 |
|
289 |
METHOD check_self; |
290 |
END check_self; |
291 |
|
292 |
METHOD check_all; |
293 |
RUN check_self; |
294 |
END check_all; |
295 |
|
296 |
METHOD bound_self; |
297 |
END bound_self; |
298 |
|
299 |
METHOD bound_all; |
300 |
RUN bound_self; |
301 |
END bound_all; |
302 |
|
303 |
METHOD scale_self; |
304 |
END scale_self; |
305 |
|
306 |
METHOD scale_all; |
307 |
RUN scale_self; |
308 |
END scale_all; |
309 |
END test_elm_kinetics; |
310 |
|
311 |
MODEL test_ord_kinetics; |
312 |
components IS_A set OF symbol_constant; |
313 |
nr IS_A set OF symbol_constant; |
314 |
Kr[nr] IS_A constant; |
315 |
active_energy[nr] IS_A constant; |
316 |
reac_T IS_A temperature; |
317 |
species[nr] IS_A set OF symbol_constant; |
318 |
nu[components][nr] IS_A integer_constant; |
319 |
order[components][nr] IS_A constant; |
320 |
conc[components] IS_A molar_density; |
321 |
|
322 |
|
323 |
components :==['methane','ethane','n_butane']; |
324 |
nr :==['reaction_1','reaction_2']; |
325 |
species['reaction_1'] :==['methane','ethane']; |
326 |
species['reaction_2'] :==['methane','ethane']; |
327 |
Kr['reaction_1'] :==2 {m^3/mol/s}; |
328 |
Kr['reaction_2'] :==3 {1/s}; |
329 |
active_energy['reaction_1'] :==70000 {J}; |
330 |
active_energy['reaction_2'] :==70000 {J}; |
331 |
|
332 |
nu['methane']['reaction_1'] :==-2; |
333 |
nu['ethane']['reaction_1'] :==1; |
334 |
nu['methane']['reaction_2'] :==2; |
335 |
nu['ethane']['reaction_2'] :==-1; |
336 |
(* define nu=0 FOR components not IN reactions*) |
337 |
FOR i IN nr CREATE |
338 |
FOR j IN components - [species[i]] CREATE |
339 |
nu[j][i] :==0; |
340 |
END FOR; |
341 |
END FOR; |
342 |
|
343 |
order['methane']['reaction_1'] :==2; |
344 |
order['ethane']['reaction_1'] :==0; |
345 |
order['methane']['reaction_2'] :==3; |
346 |
order['ethane']['reaction_2'] :==1; |
347 |
FOR i IN nr CREATE |
348 |
FOR j IN components - [species[i]] CREATE |
349 |
order[j][i] :==0; |
350 |
END FOR; |
351 |
END FOR; |
352 |
|
353 |
reac IS_A specify_kinetics(components, nr, Kr, active_energy, |
354 |
reac_T, species, nu, conc, order); |
355 |
|
356 |
METHODS |
357 |
METHOD default_self; |
358 |
RUN reac.default_self; |
359 |
END default_self; |
360 |
|
361 |
METHOD default_all; |
362 |
RUN reac.default_self; |
363 |
END default_all; |
364 |
|
365 |
METHOD check_self; |
366 |
END check_self; |
367 |
|
368 |
METHOD check_all; |
369 |
RUN check_self; |
370 |
END check_all; |
371 |
|
372 |
METHOD bound_self; |
373 |
END bound_self; |
374 |
|
375 |
METHOD bound_all; |
376 |
RUN bound_self; |
377 |
END bound_all; |
378 |
|
379 |
METHOD scale_self; |
380 |
END scale_self; |
381 |
|
382 |
METHOD scale_all; |
383 |
RUN scale_self; |
384 |
END scale_all; |
385 |
END test_ord_kinetics; |