Line data Source code
1 : ! ***********************************************************************
2 : !
3 : ! Copyright (C) 2010-2019 The MESA Team
4 : !
5 : ! This program is free software: you can redistribute it and/or modify
6 : ! it under the terms of the GNU Lesser General Public License
7 : ! as published by the Free Software Foundation,
8 : ! either version 3 of the License, or (at your option) any later version.
9 : !
10 : ! This program is distributed in the hope that it will be useful,
11 : ! but WITHOUT ANY WARRANTY; without even the implied warranty of
12 : ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.
13 : ! See the GNU Lesser General Public License for more details.
14 : !
15 : ! You should have received a copy of the GNU Lesser General Public License
16 : ! along with this program. If not, see <https://www.gnu.org/licenses/>.
17 : !
18 : ! ***********************************************************************
19 :
20 : module raw_rates
21 :
22 : use rates_def
23 : use const_def !, only: missing_value, dp
24 :
25 : implicit none
26 :
27 : abstract interface
28 : subroutine rate_fcn(tf, temp, fr, rr)
29 : use const_def, only: dp
30 : use ratelib, only: T_factors
31 : implicit none
32 : type (T_Factors) :: tf
33 : real(dp), intent(in) :: temp
34 : real(dp), intent(out) :: fr, rr
35 : end subroutine rate_fcn
36 : end interface
37 :
38 : logical, parameter :: show_rates = .false.
39 :
40 : contains
41 :
42 60024 : subroutine set_raw_rates(n, irs, temp, tf, rates, ierr)
43 : use rates_def, only : T_Factors
44 : integer, intent(in) :: n
45 : integer, intent(in) :: irs(:) ! (n) maps 1..n to reaction id
46 : real(dp), intent(in) :: temp
47 : type (T_Factors) :: tf
48 : real(dp), intent(inout) :: rates(:)
49 : integer, intent(out) :: ierr
50 : integer :: i, ir, op_err
51 : include 'formats'
52 60024 : ierr = 0
53 :
54 3665910 : do i=1,n
55 3605886 : ir = irs(i)
56 3605886 : if (ir <= 0) cycle
57 : op_err = 0
58 1785704 : call set_raw_rate(ir, temp, tf, rates(i), op_err)
59 1845728 : if (op_err /= 0) then
60 0 : ierr = op_err
61 0 : cycle
62 : end if
63 : end do
64 :
65 60024 : end subroutine set_raw_rates
66 :
67 :
68 1785706 : subroutine set_raw_rate(ir, temp, tf, raw_rate, ierr)
69 60024 : use ratelib
70 : use reaclib_eval, only: do_reaclib_lookup
71 : integer, intent(in) :: ir
72 : real(dp), intent(in) :: temp
73 : type (T_Factors) :: tf
74 : real(dp), intent(out) :: raw_rate
75 : integer, intent(out) :: ierr
76 : integer :: rir, reaclib_id_ir
77 :
78 1785706 : real(dp) :: rr
79 :
80 : include 'formats'
81 :
82 1785706 : ierr = 0
83 :
84 : ! See if the rate or its reverse is being loaded from a rate_table
85 1785706 : rir = reverse_reaction_id(ir)
86 1785706 : if(rir/=0) then
87 640749 : if (raw_rates_records(ir)% use_rate_table) then
88 : ! We want a rate from its table
89 0 : call eval_table(ir, tf, temp, raw_rate, rr, ierr)
90 0 : return
91 : end if
92 :
93 : ! We want to compute a rate from the "other" table
94 640749 : if (raw_rates_records(rir)% use_rate_table) then
95 0 : reaclib_id_ir = do_reaclib_lookup(reaction_name(ir), reaclib_rates% reaction_dict)
96 : ! if ir == a reverse rate (e.g r_o16_ga_c12) will have reaclib_id_ir == 0
97 : ! if ir == a forward rate (e.g r_c12_ag_o16) will have reaclib_id_ir /= 0
98 :
99 0 : if(reaclib_id_ir == 0 ) then
100 : ! We want a reverse rate from a forward table
101 0 : call eval_table_reverse(ir, rir, tf, temp, raw_rate, rr, ierr)
102 0 : return
103 : else
104 : ! We want a forward rate from a reverse table
105 0 : write(*,'(A)')
106 0 : write(*,*) "ERROR: Can not evaluate ",trim(reaction_name(ir)), &
107 0 : " from detailed balance with ",trim(reaction_name(rir))
108 0 : write(*,*) "Provide either both rates or only provide ",trim(reaction_name(ir))
109 0 : write(*,'(A)')
110 0 : call mesa_error(__FILE__,__LINE__)
111 0 : return
112 : end if
113 : end if
114 1144957 : else if (raw_rates_records(ir)% use_rate_table) then
115 : ! Only ir is set as a table rate and rate does not have a reverse
116 38 : call eval_table(ir, tf, temp, raw_rate, rr, ierr)
117 38 : return
118 : end if
119 :
120 :
121 10019 : select case(ir)
122 :
123 : case(ir_he4_he4_he4_to_c12) ! triple alpha to c12
124 10019 : call do1(rate_tripalf_jina)
125 :
126 : case(ir_c12_to_he4_he4_he4) ! c12 to 3 alpha
127 10019 : call do1_reverse(rate_tripalf_jina)
128 :
129 : case(ir_c12_ag_o16)
130 10019 : call do1(rate_c12ag_jina)
131 :
132 : ! case(ir_o16_ga_c12) ! o16(g, a)c12
133 : ! call do1(rate_c12ag_nacre)
134 :
135 : case(ir1212) ! c12(c12,n)mg23, c12(c12,p)na23, c12(c12,a)ne20
136 10019 : call do1(rate_c12c12_fxt_multi)
137 : ! NOTE: Gasques option for c12+c12 is implemented in net, not in rates.
138 :
139 : case(ir1216)
140 10019 : call do1(rate_c12o16_fxt)
141 :
142 : case(ir1616) ! o16 + o16 -> si28 + he4
143 10019 : call do1(rate_o16o16_jina)
144 :
145 : case(ir1216_to_mg24) ! ! c12 + o16 -> mg24 + he4
146 18 : call do1(rate_c12o16a_jina)
147 :
148 : case(ir1216_to_si28) ! ! c12 + o16 -> si28
149 18 : call do1(rate_c12o16p_jina)
150 :
151 : case(ir1616a) ! o16(o16, a)si28
152 18 : call do1(rate_o16o16a_jina)
153 :
154 : case(ir1616g) ! o16(o16, g)s32
155 : ! no jina rate
156 18 : call do1(rate_o16o16g_fxt)
157 :
158 : case(ir1616p_aux) ! o16(o16, p)p31
159 18 : call do1(rate_o16o16p_jina)
160 :
161 : case(ir1616ppa) ! o16(o16, p)p31(p, a)si28
162 18 : call do1(rate_o16o16p_jina)
163 :
164 : case(ir1616ppg) ! o16(o16, p)p31(p, g)s32
165 18 : call do1(rate_o16o16p_jina)
166 :
167 : case(ir_he3_ag_be7) ! he3(he4, g)be7
168 10019 : call do1(rate_he3he4_nacre)
169 :
170 : case(ir34_pp2) ! he4(he3, g)be7(e-, nu)li7(p, a)he4
171 10019 : call do1(rate_he3he4_nacre)
172 :
173 : case(ir34_pp3) ! he4(he3, g)be7(p, g)b8(e+, nu)be8(, a)he4
174 10019 : call do1(rate_he3he4_nacre)
175 :
176 : case(ir_b8_gp_be7) ! be7(p, g)b8
177 : ! no jina rate
178 10019 : call do1_reverse(rate_be7pg_nacre)
179 :
180 : case(ir_be7_wk_li7) ! be7(e-, nu)li7 ! This is now actually handled by a custom rate
181 10019 : call do1(rate_be7em_fxt)
182 :
183 :
184 : case(ir_c12_pg_n13)
185 10019 : call do1(rate_c12pg_nacre)
186 :
187 : case(ir_c13_pg_n14)
188 10019 : call do1(rate_c13pg_nacre)
189 :
190 : case(ir_f17_gp_o16) ! f17(g, p)o16
191 10019 : call do1_reverse(rate_o16pg_nacre)
192 :
193 : case(ir_f19_pa_o16) ! f19(p, a)o16
194 10019 : call do1(rate_f19pa_nacre)
195 :
196 : case(ir_h2_be7_to_h1_he4_he4) ! be7(d, p)2he4
197 10019 : call do1(rate_be7dp_jina)
198 :
199 : case(ir_h2_h2_to_he4) ! h2(h2, g)he4
200 10019 : call do1(rate_ddg_jina)
201 :
202 : case(ir_h2_he3_to_h1_he4) ! he3(d, p)he4
203 10019 : call do1(rate_he3d_jina)
204 :
205 : case(ir_h2_pg_he3)
206 10019 : call do1(rate_dpg_nacre)
207 :
208 : case(ir_he3_be7_to_h1_h1_he4_he4) ! be7(he3, 2p)2he4
209 10019 : call do1(rate_be7he3_jina)
210 :
211 : case(ir_he3_he3_to_h1_h1_he4) ! he3(he3, 2p)he4
212 10019 : call do1(rate_he3he3_nacre)
213 :
214 : case(ir_h1_h1_he4_to_he3_he3) ! he4(2p, he3)he3
215 18 : call do1_reverse(rate_he3he3_nacre)
216 :
217 : case(ir_li7_pa_he4) ! li7(p, a)he4
218 10019 : call do1(rate_li7pa_nacre)
219 :
220 : case(ir_he4_ap_li7) ! li7(p, a)he4
221 18 : call do1_reverse(rate_li7pa_nacre)
222 :
223 : case(ir_n13_gp_c12) ! n13(g, p)c12
224 10019 : call do1_reverse(rate_c12pg_nacre)
225 :
226 : case(ir_n13_pg_o14)
227 10019 : call do1(rate_n13pg_nacre)
228 :
229 : case(ir_n14_gp_c13) ! n14(g, p)c13
230 10019 : call do1_reverse(rate_c13pg_nacre)
231 :
232 : case(ir_n14_pg_o15)
233 10019 : call do1(rate_n14pg_nacre)
234 :
235 : case(ir_n15_ap_o18) ! n15(a, p)o18
236 10019 : call do1_reverse(rate_o18pa_nacre)
237 :
238 : case(ir_ne20_ga_o16) ! ne20(g, a)o16
239 10019 : call do1_reverse(rate_o16ag_nacre)
240 :
241 : case(ir_o14_gp_n13) ! o14(g, p)n13
242 10019 : call do1_reverse(rate_n13pg_nacre)
243 :
244 : case(ir_o15_gp_n14) ! o15(g, p)n14
245 10019 : call do1_reverse(rate_n14pg_nacre)
246 :
247 : case(ir_o16_ag_ne20) ! o16(a, g)ne20
248 10019 : call do1(rate_o16ag_nacre)
249 :
250 : case(ir_o16_ap_f19) ! o16(a, p)f19
251 10019 : call do1_reverse(rate_f19pa_nacre)
252 :
253 : case(ir_o16_pg_f17)
254 10019 : call do1(rate_o16pg_nacre)
255 :
256 : case(ir_o18_pa_n15) ! o18(p, a)n15 and n15(a, p)o18
257 10019 : call do1(rate_o18pa_nacre)
258 :
259 : case(iral27pa_aux) ! al27(p, a)mg24
260 18 : call do1_reverse(rate_mg24ap_jina)
261 :
262 : case(iral27pg_aux) ! al27(p, g)si28
263 18 : call do1(rate_al27pg_jina)
264 :
265 : case(irar36ap_aux) ! ar36(a, p)k39
266 18 : call do1(rate_ar36ap_jina)
267 :
268 : case(irar36ap_to_ca40)
269 18 : call do1(rate_ar36ap_jina)
270 :
271 : case(irar36gp_aux) ! ar36(g, p)cl35
272 18 : call do1_reverse(rate_cl35pg_jina)
273 :
274 : case(irar36gp_to_s32)
275 18 : call do1_reverse(rate_cl35pg_jina)
276 :
277 : case(irbe7ec_li7_aux) ! be7(e-, nu)li7(p, a)he4
278 10019 : call do1(rate_be7em_fxt)
279 :
280 : case(irbe7pg_b8_aux) ! be7(p, g)b8(e+, nu)be8(, a)he4
281 10019 : call do1(rate_be7pg_nacre)
282 :
283 : case(irc12_to_c13) ! c12(p, g)n13(e+nu)c13
284 18 : call do1(rate_c12pg_nacre)
285 :
286 : case(irc12_to_n14) ! c12(p, g)n13(e+nu)c13(p, g)n14
287 10019 : call do1(rate_c12pg_nacre)
288 :
289 : case(irc12ap_aux) ! c12(a, p)n15
290 10019 : call do1_reverse(rate_n15pa_jina)
291 :
292 : case(irc12ap_to_o16) ! c12(a, p)n15(p, g)o16
293 10019 : call do1_reverse(rate_n15pa_jina)
294 :
295 : case(irca40ap_aux) ! ca40(a, p)sc43
296 18 : call do1(rate_ca40ap_jina)
297 :
298 : case(irca40ap_to_ti44)
299 18 : call do1(rate_ca40ap_jina)
300 :
301 : case(irca40gp_aux) ! ca40(g, p)k39
302 18 : call do1_reverse(rate_k39pg_jina)
303 :
304 : case(irca40gp_to_ar36)
305 18 : call do1_reverse(rate_k39pg_jina)
306 :
307 : case(ircl35pa_aux) ! cl35(p, a)s32
308 18 : call do1_reverse(rate_s32ap_jina)
309 :
310 : case(ircl35pg_aux) ! cl35(p, g)ar36
311 18 : call do1(rate_cl35pg_jina)
312 :
313 : case(irco55gprot_aux) ! co55(g, prot)fe54
314 18 : call do1_reverse(rate_fe54pg_jina)
315 :
316 : case(irco55pg_aux) ! co55(p, g)ni56
317 18 : call do1(rate_co55pg_jina)
318 :
319 : case(irco55protg_aux) ! co55(prot, g)ni56
320 18 : call do1(rate_co55pg_jina)
321 :
322 : case(ircr48ap_aux) ! cr48(a, p)mn51
323 18 : call do1(rate_cr48ap_jina)
324 :
325 : case(ircr48ap_to_fe52)
326 18 : call do1(rate_cr48ap_jina)
327 :
328 : case(ircr48gp_aux) ! cr48(g, p)v47
329 18 : call do1_reverse(rate_v47pg_jina)
330 :
331 : case(ircr48gp_to_ti44)
332 18 : call do1_reverse(rate_v47pg_jina)
333 :
334 : case(irf19pg_aux) ! f19(p, g)ne20
335 10019 : call do1(rate_f19pg_jina)
336 :
337 : case(irfe52ap_aux) ! fe52(a, p)co55
338 18 : call do1(rate_fe52ap_jina)
339 :
340 : case(irfe52ap_to_ni56) ! fe52(a, p)co55(p, g)ni56
341 18 : call do1(rate_fe52ap_jina)
342 :
343 : case(irfe52aprot_aux) ! fe52(a, prot)co55
344 18 : call do1(rate_fe52ap_jina)
345 :
346 : case(irfe52aprot_to_fe54) ! fe52(a, prot)co55(g, prot)fe54
347 18 : call do1(rate_fe52ap_jina)
348 :
349 : case(irfe52aprot_to_ni56) ! fe52(a, prot)co55(prot, g)ni56
350 18 : call do1(rate_fe52ap_jina)
351 :
352 : case(irfe52gp_aux) ! fe52(g, p)mn51
353 18 : call do1_reverse(rate_mn51pg_jina)
354 :
355 : case(irfe52gp_to_cr48)
356 18 : call do1_reverse(rate_mn51pg_jina)
357 :
358 : case(irfe52neut_to_fe54) ! fe52(neut, g)fe53(neut, g)fe54
359 18 : raw_rate = -1 ! rate calculated by special routine.
360 :
361 : case(irfe52ng_aux) ! fe52(n, g)fe53
362 18 : call do1(rate_fe52ng_jina)
363 :
364 : case(irfe53gn_aux) ! fe53(g, n)fe52
365 18 : call do1_reverse(rate_fe52ng_jina)
366 :
367 : case(irfe53ng_aux) ! fe53(n, g)fe54
368 18 : call do1(rate_fe53ng_jina)
369 :
370 : case(irfe54a_to_ni56) ! fe54 + alpha -> ni56 + 2 neut
371 18 : call do1(rate_fe54a_jina)
372 :
373 : case(irfe54an_aux) ! fe54(a,n)ni57
374 18 : call do1(rate_fe54an_jina)
375 :
376 : case(irfe54an_to_ni56) ! fe54(a,n)ni57(g,n)ni56
377 18 : call do1(rate_fe54an_jina)
378 :
379 : case(irfe54aprot_to_fe56) ! fe54(a, prot)co57(g, prot)fe56
380 18 : call do1(rate_fe54ap_jina)
381 :
382 : case(irfe54g_to_fe52) ! fe54(g, neut)fe53(g, neut)fe52
383 18 : call do1_reverse(rate_fe53ng_jina)
384 :
385 : case(irfe54ng_aux) ! fe54(neut, g)fe55
386 18 : call do1(rate_fe54ng_jina)
387 :
388 : case(irfe54ng_to_fe56) ! fe54(neut, g)fe55(neut, g)fe56
389 18 : raw_rate = -1 ! rate calculated by special routine.
390 :
391 : case(irfe54prot_to_fe52) ! fe54(prot, g)co55(prot, a)fe52
392 18 : raw_rate = -1 ! rate calculated by special routine.
393 :
394 : case(irfe54prot_to_ni56) ! fe54(prot, g)co55(prot, g)ni56
395 18 : raw_rate = -1 ! rate calculated by special routine.
396 :
397 : case(irfe54protg_aux) ! fe54(prot, g)co55
398 18 : call do1(rate_fe54pg_jina)
399 :
400 : case(irfe55gn_aux) ! fe55(g, neut)fe54
401 18 : call do1_reverse(rate_fe54ng_jina)
402 :
403 : case(irfe55ng_aux) ! fe55(neut, g)fe56
404 18 : call do1(rate_fe55ng_jina)
405 :
406 : case(irfe56ec_fake_to_mn56)
407 18 : raw_rate = -1 ! rate calculated by special routine.
408 :
409 : case(irfe56ec_fake_to_mn57)
410 18 : raw_rate = -1 ! rate calculated by special routine.
411 :
412 : case(irfe56ec_fake_to_cr56)
413 18 : raw_rate = -1 ! rate calculated by special routine.
414 :
415 : case(irfe56ec_fake_to_cr57)
416 18 : raw_rate = -1 ! rate calculated by special routine.
417 :
418 : case(irfe56ec_fake_to_cr58)
419 18 : raw_rate = -1 ! rate calculated by special routine.
420 :
421 : case(irfe56ec_fake_to_cr59)
422 18 : raw_rate = -1 ! rate calculated by special routine.
423 :
424 : case(irfe56ec_fake_to_cr60)
425 18 : raw_rate = -1 ! rate calculated by special routine.
426 :
427 : case(irfe56ec_fake_to_cr61)
428 18 : raw_rate = -1 ! rate calculated by special routine.
429 :
430 : case(irfe56ec_fake_to_cr62)
431 18 : raw_rate = -1 ! rate calculated by special routine.
432 :
433 : case(irfe56ec_fake_to_cr63)
434 18 : raw_rate = -1 ! rate calculated by special routine.
435 :
436 : case(irfe56ec_fake_to_cr64)
437 18 : raw_rate = -1 ! rate calculated by special routine.
438 :
439 : case(irfe56ec_fake_to_cr65)
440 18 : raw_rate = -1 ! rate calculated by special routine.
441 :
442 : case(irfe56ec_fake_to_cr66)
443 18 : raw_rate = -1 ! rate calculated by special routine.
444 :
445 : case(irfe56ee_to_ni56)
446 18 : raw_rate = -1 ! rate calculated by special routine.
447 :
448 : case(irfe56gn_aux) ! fe56(g, neut)fe55
449 18 : call do1_reverse(rate_fe55ng_jina)
450 :
451 : case(irfe56gn_to_fe54) ! fe56(g, neut)fe55(g, neut)fe54
452 18 : raw_rate = -1 ! rate calculated by special routine.
453 :
454 : case(irfe56prot_to_fe54) ! fe56(prot, g)co57(prot, a)fe54
455 18 : raw_rate = -1 ! rate calculated by special routine.
456 :
457 : case(irh2_protg_aux) ! h2(prot, g)he3
458 18 : call do1(rate_dpg_fxt)
459 :
460 : case(irh2g_neut_aux) ! h2(g, neut)prot
461 18 : call do1_reverse(rate_png_fxt)
462 :
463 : case(irhe3_neutg_aux) ! he3(neut, g)he4
464 18 : call do1(rate_he3ng_fxt)
465 :
466 : case(irhe3gprot_aux) ! he3(g, prot)h2
467 18 : call do1_reverse(rate_dpg_fxt)
468 :
469 : case(irhe4_breakup) ! he4(g, neut)he3(g, prot)h2(g, neut)prot
470 18 : raw_rate = -1 ! rate calculated by special routine.
471 :
472 : case(irhe4_rebuild) ! prot(neut, g)h2(prot, g)he3(neut, g)he4
473 18 : raw_rate = -1 ! rate calculated by special routine.
474 :
475 : case(irhe4g_neut_aux) ! he4(g, neut)he3
476 18 : call do1_reverse(rate_he3ng_fxt) ! no jina rate
477 :
478 : case(irk39pa_aux) ! k39(p, a)ar36
479 18 : call do1_reverse(rate_ar36ap_jina)
480 :
481 : case(irk39pg_aux) ! k39(p, g)ca40
482 18 : call do1(rate_k39pg_jina)
483 :
484 : case(irmg24ap_aux) ! mg24(a, p)al27
485 18 : call do1(rate_mg24ap_jina)
486 :
487 : case(irmg24ap_to_si28)
488 18 : call do1(rate_mg24ap_jina)
489 :
490 : case(irmg24gp_aux) ! mg24(g, p)na23
491 18 : call do1_reverse(rate_na23pg_jina)
492 :
493 : case(irmg24gp_to_ne20) ! mg24(g, p)na23(p, a)ne20
494 18 : call do1_reverse(rate_na23pg_jina)
495 :
496 : case(irmn51pg_aux) ! mn51(p, g)fe52
497 18 : call do1(rate_mn51pg_jina)
498 :
499 : case(irn14_to_c12) ! n14(p, g)o15(e+nu)n15(p, a)c12
500 10019 : call do1(rate_n14pg_nacre)
501 :
502 : case(irn14_to_n15) ! n14(p, g)o15(e+nu)n15
503 18 : call do1(rate_n14pg_nacre)
504 :
505 : case(irn14_to_o16) ! n14(p, g)o15(e+nu)n15(p, g)o16
506 10019 : call do1(rate_n14pg_nacre)
507 :
508 : case(irn14ag_lite) ! n14 + 1.5 alpha => ne20
509 10019 : call do1(rate_n14ag_jina)
510 :
511 : case(irn14ag_to_o18) ! n14(a, g)f18(e+nu)o18
512 10019 : call do1(rate_n14ag_jina)
513 :
514 : case(irn14gc12) ! n14 => c12 + neut + prot
515 18 : call do1_reverse(rate_c13pg_nacre)
516 :
517 : case(irn14pg_aux) ! n14(p, g)o15
518 10019 : call do1(rate_n14pg_nacre)
519 :
520 : case(irn15pa_aux) ! n15(p, a)c12
521 10019 : call do1( rate_n15pa_jina)
522 :
523 : case(irn15pg_aux) ! n15(p, g)o16
524 10019 : call do1(rate_n15pg_jina)
525 :
526 : case(irna23pa_aux) ! na23(p, a)ne20
527 10019 : call do1(rate_na23pa_jina)
528 :
529 : case(irna23pg_aux) ! na23(p, g)mg24
530 10019 : call do1(rate_na23pg_jina)
531 :
532 : case(irne18ag_to_mg24) ! ne18(a, g)mg22 -> mg24
533 18 : call do1(rate_ne18ag_jina)
534 :
535 : case(irne18ap_to_mg22) ! ne18(a, p)na21(p, g)mg22
536 10019 : call do1(rate_ne18ap_jina)
537 :
538 : case(irne18ap_to_mg24) ! ne18(a, p)na21(p, g)mg22 -> mg24
539 18 : call do1(rate_ne18ap_jina)
540 :
541 : case(irne19pg_to_mg22) ! ne19(p, g)na20(p, g)mg21(e+nu)na21(p, g)mg22
542 10019 : call do1(rate_ne19pg_jina)
543 :
544 : case(irne19pg_to_mg24) ! ne19(p, g)na20(p, g)mg21(e+nu)na21(p, g)mg22 -> mg24
545 18 : call do1(rate_ne19pg_jina)
546 :
547 : case(irne20ap_aux) ! ne20(a, p)na23
548 10019 : call do1(rate_ne20ap_jina)
549 :
550 : case(irne20ap_to_mg24) ! ne20(a, p)na23(p, g)mg24
551 10019 : call do1(rate_ne20ap_jina)
552 :
553 : case(irne20gp_aux) ! ne20(g, p)f19
554 18 : call do1_reverse(rate_f19pg_jina)
555 :
556 : case(irne20gp_to_o16) ! ne20(g, p)f19(p, a)o16
557 18 : call do1_reverse(rate_f19pg_jina)
558 :
559 : case(irne20pg_to_mg22) ! ne20(p, g)na21(p, g)mg22
560 10019 : call do1(rate_ne20pg_nacre)
561 :
562 : case(irne20pg_to_mg24) ! ne20(p, g)na21(p, g)mg22 -> mg24
563 18 : call do1(rate_ne20pg_nacre)
564 :
565 : case(irneut_to_prot) ! neut(e+nu)prot
566 10019 : raw_rate = -1 ! rate calculated by special routine.
567 :
568 : case(irni56ec_to_fe54) ! ni56 + 2 e- => 56/54*fe54
569 18 : raw_rate = -1 ! rate calculated by special routine.
570 :
571 : case(irni56ec_to_fe56) ! ni56 + 2 e- => fe56
572 10019 : raw_rate = -1 ! rate calculated by special routine.
573 :
574 : case(irni56ec_to_co56)
575 10019 : raw_rate = -1 ! rate calculated by special routine.
576 :
577 : case(irco56ec_to_fe56)
578 10019 : raw_rate = -1 ! rate calculated by special routine.
579 :
580 : case(irni56gp_aux) ! ni56(g, p)co55
581 18 : call do1_reverse(rate_co55pg_jina)
582 :
583 : case(irni56gp_to_fe52) ! ni56(g, p)co55(p, a)fe52
584 18 : raw_rate = -1 ! rate calculated by special routine.
585 :
586 : case(irni56gprot_aux) ! ni56(g, prot)co55
587 18 : call do1_reverse(rate_co55pg_jina)
588 :
589 : case(irni56gprot_to_fe52) ! ni56(g, prot)co55(prot, a)fe52
590 18 : raw_rate = -1 ! rate calculated by special routine.
591 :
592 : case(irni56gprot_to_fe54) ! ni56(g, prot)co55(g, prot)fe54
593 18 : raw_rate = -1 ! rate calculated by special routine.
594 :
595 : case(irni56ng_to_fe54) ! ni56(n,g)ni57(n,a)fe54
596 18 : raw_rate = -1 ! rate calculated by special routine.
597 :
598 : case(irni57na_aux) ! ni57(n,a)fe54
599 18 : call do1_reverse(rate_fe54an_jina)
600 :
601 : case(iro16_to_n14) ! o16(p, g)f17(e+nu)o17(p, a)n14
602 10019 : call do1(rate_o16pg_nacre)
603 :
604 : case(iro16_to_o17) ! o16(p, g)f17(e+nu)o17
605 18 : call do1(rate_o16pg_nacre)
606 :
607 : case(iro16ap_aux) ! o16(a, p)f19
608 10019 : call do1_reverse(rate_f19pa_nacre)
609 :
610 : case(iro16ap_to_ne20) ! o16(a, p)f19(p, a)ne20
611 10019 : call do1_reverse(rate_f19pa_nacre)
612 :
613 : case(iro16gp_aux) ! o16(g, p)n15
614 10019 : call do1_reverse(rate_n15pg_jina)
615 :
616 : case(iro16gp_to_c12) ! o16(g, p)n15(p, a)c12
617 18 : call do1_reverse(rate_n15pg_jina)
618 :
619 : case(iro17_to_o18) ! o17(p, g)f18(e+nu)o18
620 18 : call do1(rate_o17pg_jina)
621 :
622 : case(irp31pa_aux) ! p31(p, a)si28
623 18 : call do1_reverse(rate_si28ap_jina)
624 :
625 : case(irp31pg_aux) ! p31(p, g)s32
626 18 : call do1(rate_p31pg_jina)
627 :
628 : case(irpep_to_he3) ! p(e-p, nu)h2(p, g)he3
629 10019 : call do1(rate_pep_fxt)
630 :
631 : case(irpp_to_he3) ! p(p, e+nu)h2(p, g)he3
632 10019 : call do1(rate_pp_nacre)
633 :
634 : case(irprot_neutg_aux) ! prot(neut, g)h2
635 18 : call do1(rate_png_fxt)
636 :
637 : case(irprot_to_neut) ! prot(e-nu)neut
638 10019 : raw_rate = -1 ! rate calculated by special routine.
639 :
640 : case(irs32ap_aux) ! s32(a, p)cl35
641 18 : call do1(rate_s32ap_jina)
642 :
643 : case(irs32ap_to_ar36)
644 18 : call do1(rate_s32ap_jina)
645 :
646 : case(irs32gp_aux) ! s32(g, p)p31
647 18 : call do1_reverse(rate_p31pg_jina)
648 :
649 : case(irs32gp_to_si28)
650 18 : call do1_reverse(rate_p31pg_jina)
651 :
652 : case(irsc43pa_aux) ! sc43(p, a)ca40
653 18 : call do1_reverse(rate_ca40ap_jina)
654 :
655 : case(irsc43pg_aux) ! sc43(p, g)ti44
656 18 : call do1(rate_sc43pg_jina)
657 :
658 : case(irsi28ap_aux) ! si28(a, p)p31
659 18 : call do1(rate_si28ap_jina)
660 :
661 : case(irsi28ap_to_s32)
662 18 : call do1(rate_si28ap_jina)
663 :
664 : case(irsi28gp_aux) ! si28(g, p)al27
665 18 : call do1_reverse(rate_al27pg_jina)
666 :
667 : case(irsi28gp_to_mg24)
668 18 : call do1_reverse(rate_al27pg_jina)
669 :
670 : case(irti44ap_aux) ! ti44(a, p)v47
671 18 : call do1(rate_ti44ap_jina)
672 :
673 : case(irti44ap_to_cr48)
674 18 : call do1(rate_ti44ap_jina)
675 :
676 : case(irti44gp_aux) ! ti44(g, p)sc43
677 18 : call do1_reverse(rate_sc43pg_jina)
678 :
679 : case(irti44gp_to_ca40)
680 18 : call do1_reverse(rate_sc43pg_jina)
681 :
682 : case(irv47pa_aux) ! v47(p, a)ti44
683 18 : call do1_reverse(rate_ti44ap_jina)
684 :
685 : case(irv47pg_aux) ! v47(p, g)cr48
686 18 : call do1(rate_v47pg_jina)
687 :
688 : case(ir_h1_h1_wk_h2) ! p(p, e+nu)h2
689 10019 : call do1(rate_pp_nacre)
690 :
691 : case(ir_h1_h1_ec_h2) ! p(e-p, nu)h2
692 10019 : call do1(rate_pep_fxt)
693 :
694 : case(irn14ag_to_ne22) ! n14(a, g)f18(e+nu)o18(a, g)ne22
695 18 : call do1(rate_n14ag_jina)
696 :
697 : case(irf19pa_aux) ! f19(p, a)o16
698 10019 : call do1(rate_f19pa_nacre)
699 :
700 : case(ir_be7_pg_b8)
701 10019 : call do1(rate_be7pg_nacre)
702 :
703 : case(ir_b8_wk_he4_he4) ! b8(p=>n)be8=>2 he4
704 10019 : call do1(rate_b8ep)
705 :
706 : case(irmn51pa_aux) ! mn51(p, a)cr48
707 18 : call do1_reverse(rate_cr48ap_jina)
708 :
709 : case(irfe54gn_aux) ! fe54(g, n)fe53
710 18 : call do1_reverse(rate_fe53ng_jina)
711 :
712 : case(irco55pa_aux) ! co55(p, a)fe52
713 18 : call do1_reverse(rate_fe52ap_jina)
714 :
715 : case(irco55prota_aux) ! co55(prot, a)fe52
716 18 : call do1_reverse(rate_fe52ap_jina)
717 :
718 : case(ir_h1_he3_wk_he4) ! he3(p, e+nu)he4
719 10019 : call do1(rate_hep_fxt)
720 :
721 : case(ir_he3_ng_he4)
722 10019 : call do1(rate_he3ng_fxt)
723 :
724 : case(ir_he4_gn_he3)
725 10019 : call do1_reverse(rate_he3ng_fxt)
726 :
727 : case(ir_h1_ng_h2)
728 10019 : call do1(rate_png_fxt)
729 :
730 : case(ir_h2_gn_h1)
731 10019 : call do1_reverse(rate_png_fxt)
732 :
733 : case(ir_he3_gp_h2)
734 10019 : call do1_reverse(rate_dpg_nacre)
735 :
736 : case(ir_c12_c12_to_h1_na23)
737 18 : call do1(rate_c12_c12_to_h1_na23_jina)
738 :
739 : case(ir_he4_ne20_to_c12_c12)
740 18 : call do1(rate_he4_ne20_to_c12_c12_jina)
741 :
742 : case(ir_c12_c12_to_he4_ne20)
743 18 : call do1_reverse(rate_he4_ne20_to_c12_c12_jina)
744 :
745 : case(ir_he4_mg24_to_c12_o16)
746 18 : call do1(rate_he4_mg24_to_c12_o16_jina)
747 :
748 : case default
749 1785668 : call do_default(ierr)
750 :
751 : end select
752 :
753 :
754 1785706 : if(associated(rates_other_rate_get)) then
755 0 : call rates_other_rate_get(ir, temp, tf, raw_rate, ierr)
756 : end if
757 :
758 :
759 : contains
760 :
761 :
762 1021884 : subroutine do_default(ierr)
763 : integer, intent(out) :: ierr ! set ierr to -1 if cannot find rate
764 : real(dp) :: lambda, dlambda_dlnT, rlambda, drlambda_dlnT
765 : include 'formats'
766 : ierr = 0
767 : ! look for rate in reaclib
768 : call get_reaclib_rate_and_dlnT( &
769 1021884 : ir, temp, lambda, dlambda_dlnT, rlambda, drlambda_dlnT, ierr)
770 1021884 : raw_rate = lambda
771 1785706 : end subroutine do_default
772 :
773 :
774 532069 : subroutine do1(rate_fcn1)
775 : procedure(rate_fcn) :: rate_fcn1
776 532069 : call eval_raw_rate(ir, rate_fcn1, tf, temp, raw_rate, rr, ierr)
777 532069 : end subroutine do1
778 :
779 :
780 181134 : subroutine do1_reverse(rate_fcn1)
781 : procedure(rate_fcn) :: rate_fcn1
782 181134 : call eval_raw_rate(ir, rate_fcn1, tf, temp, rr, raw_rate, ierr)
783 181134 : end subroutine do1_reverse
784 :
785 : end subroutine set_raw_rate
786 :
787 :
788 0 : subroutine rate_fcn_null(tf, temp, fr, rr)
789 : use ratelib, only: T_factors
790 : type (T_Factors) :: tf
791 : real(dp), intent(in) :: temp
792 : real(dp), intent(out) :: fr, rr
793 0 : fr = -1; rr = -1
794 0 : end subroutine rate_fcn_null
795 :
796 :
797 0 : subroutine eval_which_raw_rate( &
798 : ir, which_rate, &
799 : rate_fcn1, rate_fcn2, rate_fcn3, rate_fcn4, &
800 : tf, temp, fr, rr, ierr)
801 0 : use interp_1d_lib, only: interp_values
802 : use ratelib
803 : use const_def, only: pi
804 : use math_lib
805 : integer, intent(in) :: ir, which_rate
806 : procedure(rate_fcn) :: rate_fcn1, rate_fcn2, rate_fcn3, rate_fcn4
807 : type (T_Factors) :: tf
808 : real(dp), intent(in) :: temp
809 : real(dp), intent(out) :: fr, rr
810 : integer, intent(out) :: ierr
811 : include 'formats'
812 : ierr = 0
813 :
814 0 : call eval_raw_rate(ir, rate_fcn1, tf, temp, fr, rr, ierr)
815 :
816 0 : end subroutine eval_which_raw_rate
817 :
818 :
819 713203 : subroutine eval_raw_rate(ir, rate_fcn1, tf, temp, fr, rr, ierr)
820 0 : use ratelib
821 : integer, intent(in) :: ir ! reaction id
822 : procedure(rate_fcn) :: rate_fcn1
823 : type (T_Factors) :: tf
824 : real(dp), intent(in) :: temp
825 : real(dp), intent(out) :: fr, rr
826 : integer, intent(out) :: ierr
827 : include 'formats'
828 713203 : ierr = 0
829 713203 : call rate_fcn1(tf, temp, fr, rr)
830 713203 : if (fr < 0 .or. rr < 0) then
831 0 : write(*,1) 'invalid which rate for ' // trim(reaction_Name(ir)), fr, rr
832 0 : ierr = -1
833 : return
834 : end if
835 : if (fr == missing_value .or. rr == missing_value) then
836 : write(*,1) 'missing value for ' // trim(reaction_Name(ir)), fr, rr
837 : ierr = -1
838 : return
839 : end if
840 713203 : end subroutine eval_raw_rate
841 :
842 38 : subroutine eval_table(ir, tf, temp, fr, rr, ierr)
843 713203 : use interp_1d_lib, only: interp_values
844 : use ratelib
845 : integer, intent(in) :: ir ! reaction id
846 : type (T_Factors) :: tf
847 : real(dp), intent(in) :: temp
848 : real(dp), intent(out) :: fr, rr
849 : integer, intent(out) :: ierr
850 : integer, parameter :: nv = 1
851 114 : real(dp) :: x(nv), vals(nv)
852 : type (rate_table_info), pointer :: ri=> null()
853 : include 'formats'
854 38 : ierr = 0
855 :
856 38 : ri => raw_rates_records(ir)
857 38 : if (.not. ri% use_rate_table) then
858 0 : ierr = -1
859 0 : return
860 : end if
861 38 : if (ri% need_to_read) then
862 12 : !$omp critical (load_rate_table)
863 6 : if (ri% need_to_read) then
864 6 : call get_interp_table(ri% rate_fname, ri% nT8s, ri% T8s, ri% f1, ierr)
865 6 : ri% need_to_read = .false.
866 : end if
867 : !$omp end critical (load_rate_table)
868 6 : if (ierr /= 0) then
869 0 : write(*,*) 'failed to load table ' // trim(ri% rate_fname)
870 0 : return
871 : end if
872 : end if
873 38 : x(1) = temp*1d-8
874 38 : call interp_values(ri% T8s, ri% nT8s, ri% f1, nv, x, vals, ierr)
875 38 : fr = vals(1)
876 38 : rr = 0 ! no reverse rates for tables
877 38 : end subroutine eval_table
878 :
879 :
880 0 : subroutine eval_table_reverse(ir, rir, tf, temp, fr, rr, ierr)
881 38 : use interp_1d_lib, only: interp_values
882 : use ratelib
883 : use reaclib_eval, only: compute_some_inverse_lambdas,&
884 : do_reaclib_indices_for_reaction
885 : integer, intent(in) :: ir, rir ! reaction id, reverse reaction id
886 : type (T_Factors) :: tf
887 : real(dp), intent(in) :: temp
888 : real(dp), intent(out) :: fr, rr
889 : integer, intent(out) :: ierr
890 0 : real(dp), dimension(1):: ln_lambda, lambda, dlambda_dlnT, &
891 0 : inv_lambda, dinv_lambda_dlnT
892 : integer :: lo, hi
893 : real(dp) :: fr_table
894 : include 'formats'
895 : ierr = 0
896 :
897 0 : call eval_table(rir, tf, temp, fr_table, rr, ierr)
898 :
899 0 : lambda = 1d0
900 0 : ln_lambda = 0d0
901 0 : dlambda_dlnT = 0d0
902 :
903 0 : call do_reaclib_indices_for_reaction(reaction_name(rir), reaclib_rates, lo, hi, ierr)
904 0 : if(ierr/=0) then
905 0 : write(*,*) "Error: Could not get reaclib index for rate ",rir, ' ',trim(reaction_name(rir))
906 : return
907 : end if
908 :
909 0 : inv_lambda = 0d0
910 0 : dinv_lambda_dlnT = 0d0
911 :
912 : call compute_some_inverse_lambdas(1, lo, lo, &
913 : tf%T9, reaclib_rates, &
914 : ln_lambda, lambda, dlambda_dlnT, &
915 0 : inv_lambda, dinv_lambda_dlnT)
916 :
917 0 : fr = inv_lambda(1) * fr_table
918 :
919 0 : rr = 0
920 0 : end subroutine eval_table_reverse
921 :
922 :
923 6 : subroutine get_interp_table(f_name, nT8s, T8s_out, f1_out, ierr)
924 0 : use interp_1d_def, only: pm_work_size
925 : use interp_1d_lib, only: interp_pm
926 : use utils_lib
927 : use math_lib, only: str_to_vector
928 : character (len=*), intent(in) :: f_name
929 : integer, intent(out) :: nT8s
930 : real(dp), pointer :: T8s_out(:) ! will be allocated. (nT8s)
931 : real(dp), pointer :: f1_out(:) ! will be allocated. (4,nT8s)
932 : integer, intent(out) :: ierr
933 :
934 : integer :: iounit, j, nvec
935 : real(dp), pointer :: work(:)=> null()
936 : real(dp), pointer :: T8s(:)=> null()
937 : real(dp), pointer :: f1(:)=> null(), f(:,:)=> null()
938 : character (len=256) :: line, rate_file
939 126 : real(dp), target :: vec_ary(20)
940 : real(dp), pointer :: vec(:)=> null()
941 :
942 6 : ierr = 0
943 6 : vec => vec_ary
944 :
945 :
946 : ! Look for the file based on its name first
947 :
948 6 : rate_file = trim(f_name)
949 :
950 6 : open(newunit=iounit,file=trim(rate_file),action='read',status='old',iostat=ierr)
951 6 : if (ierr /= 0) then
952 : ! Look in rates_table_dir
953 6 : rate_file = trim(rates_table_dir) // '/' // trim(f_name)
954 6 : open(newunit=iounit,file=trim(rate_file),action='read',status='old',iostat=ierr)
955 6 : if (ierr /= 0) then
956 0 : write(*,*) 'ERROR: cannot open rate info file ' // trim(rate_file)
957 : !return
958 0 : call mesa_error(__FILE__,__LINE__)
959 : end if
960 : end if
961 :
962 : do ! read until reach line starting with an integer (nT8s)
963 28 : ierr = 0
964 28 : read(iounit, fmt=*, iostat=ierr) nT8s
965 28 : if (ierr == 0) exit
966 28 : if (is_iostat_end(ierr)) exit
967 : end do
968 6 : if (failed('Failed to find line starting with an integer for nT8s')) return
969 :
970 6 : allocate(T8s(nT8s), f1(4*nT8s), stat=ierr)
971 6 : if (failed('allocate')) return
972 6 : f(1:4,1:nT8s) => f1(1:4*nT8s)
973 :
974 770 : do j=1,nT8s
975 764 : read(iounit,'(a)',iostat=ierr) line
976 764 : if (ierr == 0) call str_to_vector(line, vec, nvec, ierr)
977 764 : if (failed('read table')) return
978 764 : if (nvec < 2) then
979 0 : ierr = -1
980 0 : if (failed('read table')) return
981 : end if
982 764 : T8s(j) = vec(1)
983 770 : f(1,j) = vec(2)
984 : end do
985 :
986 6 : allocate(work(nT8s*pm_work_size), stat=ierr)
987 6 : if (failed('allocate')) return
988 :
989 : call interp_pm(T8s, nT8s, f1, pm_work_size, work, &
990 6 : 'rates get_interp_table', ierr)
991 6 : deallocate(work)
992 :
993 6 : if (failed('interp_pm')) return
994 :
995 6 : close(iounit)
996 :
997 : ! don't set the pointers until have finished setting up the data
998 :
999 6 : if (associated(T8s_out)) deallocate(T8s_out)
1000 6 : if (associated(f1_out)) deallocate(f1_out)
1001 :
1002 6 : T8s_out => T8s
1003 12 : f1_out => f1
1004 :
1005 : contains
1006 :
1007 788 : logical function failed(str)
1008 : character (len=*), intent(in) :: str
1009 788 : failed = .false.
1010 788 : if (ierr == 0) return
1011 0 : close(iounit)
1012 0 : write(*,*) trim(str) // ' failed in reading ' // trim(rate_file)
1013 0 : failed = .true.
1014 0 : return
1015 6 : end function failed
1016 :
1017 : end subroutine get_interp_table
1018 :
1019 :
1020 1021884 : subroutine get_reaclib_rate_and_dlnT( &
1021 : ir, temp, lambda, dlambda_dlnT, rlambda, drlambda_dlnT, ierr)
1022 : use ratelib
1023 : integer, intent(in) :: ir
1024 : real(dp), intent(in) :: temp
1025 : real(dp), intent(out) :: lambda, dlambda_dlnT, rlambda, drlambda_dlnT
1026 : integer, intent(out) :: ierr
1027 : integer :: reverse
1028 :
1029 : include 'formats'
1030 1021884 : ierr = 0
1031 :
1032 1021884 : reverse = reaction_is_reverse(ir)
1033 1021884 : if (reverse == 0) then ! that means don't know
1034 300 : reverse = reaclib_reverse(reaction_Name(ir))
1035 300 : if (reverse > 0) then
1036 135 : reaction_is_reverse(ir) = reverse
1037 : else
1038 165 : reaction_is_reverse(ir) = -1
1039 : end if
1040 1021584 : else if (reverse == -1) then ! that means reaclib_reverse returned 0
1041 : reverse = 0
1042 : end if
1043 451020 : if (reverse > 0) then
1044 : !write(*,1) 'call do_jina_reaclib_reverse ' // trim(reaction_Name(ir))
1045 450855 : call do_jina_reaclib_reverse(reaclib_rates% reaction_handle(reverse))
1046 : else
1047 : !write(*,1) 'call do_jina_reaclib ' // trim(reaction_Name(ir))
1048 571029 : call do_jina_reaclib
1049 : end if
1050 :
1051 : return
1052 :
1053 : write(*,1) 'temp', temp
1054 : write(*,1) 'lambda', lambda
1055 : write(*,1) 'dlambda_dlnT', dlambda_dlnT
1056 : write(*,1) 'rlambda', rlambda
1057 : write(*,1) 'drlambda_dlnT', drlambda_dlnT
1058 :
1059 1021884 : call mesa_error(__FILE__,__LINE__,'get_reaclib_rate_and_dlnT')
1060 :
1061 : contains
1062 :
1063 :
1064 1021884 : subroutine get_reaclib_lo_hi(ir, handle, lo, hi, ierr)
1065 1021884 : use reaclib_eval, only: do_reaclib_indices_for_reaction
1066 : integer, intent(in) :: ir
1067 : character (len=*) :: handle
1068 : integer, intent(out) :: lo, hi, ierr
1069 : include 'formats'
1070 1021884 : ierr = 0
1071 1021884 : lo = reaction_reaclib_lo(ir)
1072 1021884 : hi = reaction_reaclib_hi(ir)
1073 1021884 : if (lo > 0 .and. hi > 0) return
1074 : call do_reaclib_indices_for_reaction( &
1075 168 : handle, reaclib_rates, lo, hi, ierr)
1076 168 : if (ierr /= 0) return
1077 168 : if (lo <= 0 .or. hi <= 0) ierr = -1
1078 168 : reaction_reaclib_lo(ir) = lo
1079 168 : reaction_reaclib_hi(ir) = hi
1080 1021884 : end subroutine get_reaclib_lo_hi
1081 :
1082 :
1083 571029 : subroutine do_jina_reaclib
1084 : integer :: ierr, lo, hi
1085 : include 'formats'
1086 : ierr = 0
1087 571029 : call get_reaclib_lo_hi(ir, reaction_Name(ir), lo, hi, ierr)
1088 571029 : if (ierr /= 0) then
1089 : write(*,'(a,3x,i5)') &
1090 0 : trim(reaction_Name(ir)) // ' failed in do_jina_reaclib', ir
1091 : !call mesa_error(__FILE__,__LINE__,'raw_rates')
1092 0 : return
1093 : end if
1094 : !write(*,3) trim(reaction_Name(ir)) // ' lo hi', lo, hi
1095 : call reaclib_rate_and_dlnT( &
1096 : lo, hi, reaction_Name(ir), temp*1d-9, &
1097 571029 : lambda, dlambda_dlnT, rlambda, drlambda_dlnT, ierr)
1098 : !write(*,2) trim(reaction_Name(ir)) // ' lambda', ir, lambda
1099 571029 : if (ierr /= 0) then
1100 : write(*,'(a,3x,i5)') &
1101 0 : trim(reaction_Name(ir)) // ' failed in get_reaclib_rate_and_dlnT', ir
1102 0 : return
1103 : end if
1104 1021884 : end subroutine do_jina_reaclib
1105 :
1106 :
1107 901710 : subroutine do_jina_reaclib_reverse(reverse_handle)
1108 : character (len=*) :: reverse_handle
1109 : integer :: ierr, lo, hi, r_id
1110 : include 'formats'
1111 450855 : ierr = 0
1112 450855 : r_id = reverse_reaction_id(ir)
1113 450855 : if (r_id == 0) then ! don't know
1114 125 : r_id = get_rates_reaction_id(reverse_handle)
1115 125 : if (r_id == 0) then
1116 : write(*,'(a,3x,i5)') &
1117 0 : trim(reverse_handle) // ' failed in reaclib_index', r_id
1118 : !call mesa_error(__FILE__,__LINE__,'raw_rates')
1119 : end if
1120 125 : reverse_reaction_id(ir) = r_id
1121 : end if
1122 450855 : call get_reaclib_lo_hi(r_id, reverse_handle, lo, hi, ierr)
1123 450855 : if (ierr /= 0) then
1124 : write(*,'(a,3x,i5)') &
1125 0 : trim(reverse_handle) // ' failed in do_jina_reaclib_reverse', r_id
1126 0 : call mesa_error(__FILE__,__LINE__,'raw_rates')
1127 0 : return
1128 : end if
1129 : call reaclib_rate_and_dlnT( &
1130 : lo, hi, reverse_handle, temp*1d-9, &
1131 450855 : rlambda, drlambda_dlnT, lambda, dlambda_dlnT, ierr)
1132 450855 : if (ierr /= 0) then
1133 0 : write(*,'(a)') trim(reverse_handle) // ' failed in get_reaclib_rate_and_dlnT'
1134 0 : return
1135 : end if
1136 : end subroutine do_jina_reaclib_reverse
1137 :
1138 : end subroutine get_reaclib_rate_and_dlnT
1139 :
1140 : end module raw_rates
|