Line data Source code
1 : ! ***********************************************************************
2 : !
3 : ! Copyright (C) 2011-2019 Aaron Dotter, Ed Brown & 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 chem_def
21 :
22 : use utils_def, only: integer_dict
23 : use const_def, only: dp
24 : use math_lib, only: exp10
25 : use utils_lib, only: mesa_error
26 :
27 : implicit none
28 :
29 : ! Some notes on solar abundance scales:
30 : !
31 : ! published tabulations typically include both photospheric (spectroscopic)
32 : ! and meteoritic abundance measurements. many elements can be measured
33 : ! in both. some elements can only be measured one way (or the other). The
34 : ! following table, compiled by Frank Timmes, gives the breakdown of the
35 : ! elemental abundances for each of the solar abundance patterns included in
36 : ! MESA.
37 : !
38 : ! m = meteoritic, generally chondrites
39 : ! s = spectroscopic photosphere
40 : ! a = average of m and s
41 : ! r = recommended
42 : !
43 : !name what it is what mesa adopts
44 : !AG89 m,s m, with some s
45 : !GS98 m,s m, with some s
46 : !L03 m,s,a as tabulated
47 : !AGS05 m,s s, with some m
48 : !AGSS09 m,s s, with some m
49 : !L09 m,s,r r
50 : !A09_Prz modified AGSS09 - he, c, n, o, ne, mg, al, si, s, ar, fe.
51 : !MB22 m,s s, with some m, supplemented with agss09
52 : !AAG21 m,s s, with some m
53 :
54 : ! storage for Lodders (2003) isotopic percentages
55 : type isotopic_abundance_table_type
56 : type (integer_dict), pointer :: name_dict
57 : real(dp), dimension(:), allocatable :: isotopic_percent
58 : end type isotopic_abundance_table_type
59 :
60 : type(isotopic_abundance_table_type) :: lodders03_tab6
61 :
62 :
63 : integer, parameter :: iso_name_length = 8 ! no. characters in an nuclide name
64 : ! 2 for element, 3 for A, and 3 for isomeric state as '_' + 1 or 2 digit integer 0 to 99)
65 : integer, parameter :: long_name_length = 16 ! no. characters in a long name
66 : integer, parameter :: npart = 24 ! no. entries in partition fcn table
67 : integer, parameter :: max_el_z = 112 ! max. Z in the winvn database
68 : integer, parameter :: nhname = 3 ! no. isotopes for h in the winvn database
69 :
70 : logical, parameter :: convert_mass_excess_to_binding_energy = .true.
71 : integer, parameter :: nuclide_not_found = -1 ! warning flag
72 :
73 : integer, dimension(0:max_el_z) :: element_min_N, element_max_N
74 : ! for isos included in chem_isos
75 :
76 : ! element names
77 : character(len=iso_name_length), dimension(0:max_el_z) :: &
78 : el_name = [character(len=iso_name_length) :: &
79 : 'neut','h','he','li','be','b','c','n','o','f','ne', &
80 : 'na','mg','al','si','p','s','cl','ar','k','ca', &
81 : 'sc','ti','v','cr','mn','fe','co','ni','cu','zn', &
82 : 'ga','ge','as','se','br','kr','rb','sr','y','zr', &
83 : 'nb','mo','tc','ru','rh','pd','ag','cd','in','sn', &
84 : 'sb','te','i','xe','cs','ba','la','ce','pr','nd', &
85 : 'pm','sm','eu','gd','tb','dy','ho','er','tm','yb', &
86 : 'lu','hf','ta','w','re','os','ir','pt','au','hg', &
87 : 'tl','pb','bi','po','at','rn','fr','ra','ac','th', &
88 : 'pa','u','np','pu','am','cm','bk','cf','es','fm','md', &
89 : 'no','lr','rf','db','sg','bh','hs','mt','ds','rg','cn']
90 :
91 : character(len=long_name_length), dimension(0:max_el_z) :: &
92 : el_long_name = [character(len=long_name_length) :: &
93 : 'neutron','hydrogen','helium','lithium','beryllium', &
94 : 'boron','carbon','nitrogen','oxygen','fluorine','neon', &
95 : 'sodium','magnesium','aluminum','silicon','phosphorus', &
96 : 'sulfur','chlorine','argon','potassium','calcium', &
97 : 'scandium','titanium','vanadium','chromium','manganese', &
98 : 'iron','cobalt','nickel','copper','zinc','gallium', &
99 : 'germanium','arsenic','selenium','bromine','krypton', &
100 : 'rubidium','strontium','yttrium','zirconium','niobium', &
101 : 'molybdenum','technetium','ruthenium','rhodium', &
102 : 'palladium','silver','cadmium','indium','tin','antimony', &
103 : 'tellurium','iodine','xenon','cesium','barium', &
104 : 'lanthanum','cerium','praseodymium','neodymium','promethium', &
105 : 'samarium','europium','gadolinium','terbium','dysprosium', &
106 : 'holmium','erbium','thulium','ytterbium','lutetium', &
107 : 'hafnium','tantalum','tungsten','rhenium','osmium', &
108 : 'iridium','platinum','gold','mercury','thallium','lead', &
109 : 'bisumth','polonium','astatine','radon','francium', &
110 : 'radium','actinium','thorium','protactinium','uranium', &
111 : 'neptunium','plutonium','americium','curium','berkelium', &
112 : 'californium','einsteinium','fermium','mendelevium', &
113 : 'nobelium','lawrencium','rutherfordium','dubnium', &
114 : 'seaborgium','bohrium','hassium','meitnerium','darmstadtium', &
115 : 'roentgenium','copernicium' ]
116 :
117 : ! aluminum isomers
118 : character(len=iso_name_length), dimension(2:3) :: &
119 : al_isomers = [character(len=iso_name_length) ::'al-6','al*6']
120 : character(len=long_name_length), dimension(2:3) :: &
121 : long_al_isomers = [character(len=long_name_length) :: &
122 : 'Aluminum-gs','Aluminum-ex']
123 :
124 :
125 : ! chem element id numbers (up to Cn)
126 : ! note: for isotope i, the element id number = chem_Z(i)
127 :
128 : !periodic table, row 1
129 : integer, parameter :: e_h = 1
130 : integer, parameter :: e_he = 2
131 :
132 : !periodic table, row 2
133 : integer, parameter :: e_li = 3
134 : integer, parameter :: e_be = 4
135 : integer, parameter :: e_b = 5
136 : integer, parameter :: e_c = 6
137 : integer, parameter :: e_n = 7
138 : integer, parameter :: e_o = 8
139 : integer, parameter :: e_f = 9
140 : integer, parameter :: e_ne = 10
141 :
142 : !periodic table, row 3
143 : integer, parameter :: e_na = 11
144 : integer, parameter :: e_mg = 12
145 : integer, parameter :: e_al = 13
146 : integer, parameter :: e_si = 14
147 : integer, parameter :: e_p = 15
148 : integer, parameter :: e_s = 16
149 : integer, parameter :: e_cl = 17
150 : integer, parameter :: e_ar = 18
151 :
152 : !periodic table, row 4
153 : integer, parameter :: e_k = 19
154 : integer, parameter :: e_ca = 20
155 : integer, parameter :: e_sc = 21
156 : integer, parameter :: e_ti = 22
157 : integer, parameter :: e_v = 23
158 : integer, parameter :: e_cr = 24
159 : integer, parameter :: e_mn = 25
160 : integer, parameter :: e_fe = 26
161 : integer, parameter :: e_co = 27
162 : integer, parameter :: e_ni = 28
163 : integer, parameter :: e_cu = 29
164 : integer, parameter :: e_zn = 30
165 : integer, parameter :: e_ga = 31
166 : integer, parameter :: e_ge = 32
167 : integer, parameter :: e_as = 33
168 : integer, parameter :: e_se = 34
169 : integer, parameter :: e_br = 35
170 : integer, parameter :: e_kr = 36
171 :
172 : !periodic table, row 5
173 : integer, parameter :: e_rb = 37
174 : integer, parameter :: e_sr = 38
175 : integer, parameter :: e_y = 39
176 : integer, parameter :: e_zr = 40
177 : integer, parameter :: e_nb = 41
178 : integer, parameter :: e_mo = 42
179 : integer, parameter :: e_tc = 43
180 : integer, parameter :: e_ru = 44
181 : integer, parameter :: e_rh = 45
182 : integer, parameter :: e_pd = 46
183 : integer, parameter :: e_ag = 47
184 : integer, parameter :: e_cd = 48
185 : integer, parameter :: e_in = 49
186 : integer, parameter :: e_sn = 50
187 : integer, parameter :: e_sb = 51
188 : integer, parameter :: e_te = 52
189 : integer, parameter :: e_i = 53
190 : integer, parameter :: e_xe = 54
191 :
192 : !periodic table, row 6
193 : integer, parameter :: e_cs = 55
194 : integer, parameter :: e_ba = 56
195 : integer, parameter :: e_la = 57
196 : integer, parameter :: e_ce = 58
197 : integer, parameter :: e_pr = 59
198 : integer, parameter :: e_nd = 60
199 : integer, parameter :: e_pm = 61
200 : integer, parameter :: e_sm = 62
201 : integer, parameter :: e_eu = 63
202 : integer, parameter :: e_gd = 64
203 : integer, parameter :: e_tb = 65
204 : integer, parameter :: e_dy = 66
205 : integer, parameter :: e_ho = 67
206 : integer, parameter :: e_er = 68
207 : integer, parameter :: e_tm = 69
208 : integer, parameter :: e_yb = 70
209 : integer, parameter :: e_lu = 71
210 : integer, parameter :: e_hf = 72
211 : integer, parameter :: e_ta = 73
212 : integer, parameter :: e_w = 74
213 : integer, parameter :: e_re = 75
214 : integer, parameter :: e_os = 76
215 : integer, parameter :: e_ir = 77
216 : integer, parameter :: e_pt = 78
217 : integer, parameter :: e_au = 79
218 : integer, parameter :: e_hg = 80
219 : integer, parameter :: e_tl = 81
220 : integer, parameter :: e_pb = 82
221 : integer, parameter :: e_bi = 83
222 : integer, parameter :: e_po = 84
223 : integer, parameter :: e_at = 85
224 : integer, parameter :: e_rn = 86
225 :
226 : !periodic table, row 7
227 : integer, parameter :: e_fr = 87
228 : integer, parameter :: e_ra = 88
229 : integer, parameter :: e_ac = 89
230 : integer, parameter :: e_th = 90
231 : integer, parameter :: e_pa = 91
232 : integer, parameter :: e_u = 92
233 : integer, parameter :: e_np = 93
234 : integer, parameter :: e_pu = 94
235 : integer, parameter :: e_am = 95
236 : integer, parameter :: e_cm = 96
237 : integer, parameter :: e_bk = 97
238 : integer, parameter :: e_cf = 98
239 : integer, parameter :: e_es = 99
240 : integer, parameter :: e_fm = 100
241 : integer, parameter :: e_md = 101
242 : integer, parameter :: e_no = 102
243 : integer, parameter :: e_lr = 103
244 : integer, parameter :: e_rf = 104
245 : integer, parameter :: e_db = 105
246 : integer, parameter :: e_sg = 106
247 : integer, parameter :: e_bh = 107
248 : integer, parameter :: e_hs = 108
249 : integer, parameter :: e_mt = 109
250 : integer, parameter :: e_ds = 110
251 : integer, parameter :: e_rg = 111
252 : integer, parameter :: e_cn = 112
253 :
254 : integer, parameter :: num_chem_elements = max_el_z
255 :
256 :
257 : ! anders & grevesse 1989
258 : integer, parameter :: solsiz = 286
259 : integer, parameter :: solnamelen = 5
260 : character (len=solnamelen) :: namsol(solsiz)
261 : integer :: izsol(solsiz),iasol(solsiz),jcode(solsiz)
262 : real(dp) :: solx(solsiz), zsol, yesol ! according to AG89
263 : type (integer_dict), pointer :: Xsol_names_dict
264 :
265 :
266 : ! various values for current solar Z and Y (at photosphere)
267 : ! note that these have been reduced by diffusion from pre-MS values.
268 : ! values updated from Asplund et al. ARAA, 2009, 47, 481
269 :
270 : real(dp), parameter :: AG89_zsol = 0.0201d0
271 : real(dp), parameter :: GN93_zsol = 0.0179d0
272 : real(dp), parameter :: GS98_zsol = 0.0169d0
273 : real(dp), parameter :: L03_zsol = 0.0133d0
274 : real(dp), parameter :: AGS05_zsol = 0.0122d0
275 : real(dp), parameter :: AGSS09_zsol = 0.0134d0
276 : real(dp), parameter :: L09_zsol = 0.0154d0
277 : real(dp), parameter :: A09_Prz_zsol = 0.014d0
278 : real(dp), parameter :: MB22_photospheric_zsol = 0.0169d0
279 : real(dp), parameter :: AAG21_photospheric_zsol = 0.0139d0
280 :
281 : real(dp), parameter :: AG89_ysol = 0.2485d0
282 : real(dp), parameter :: GN93_ysol = 0.2485d0
283 : real(dp), parameter :: GS98_ysol = 0.2485d0
284 : real(dp), parameter :: L03_ysol = 0.2377d0
285 : real(dp), parameter :: AGS05_ysol = 0.2485d0
286 : real(dp), parameter :: AGSS09_ysol = 0.2485d0
287 : real(dp), parameter :: L09_ysol = 0.2751d0
288 : real(dp), parameter :: A09_Prz_ysol = 0.276d0
289 : real(dp), parameter :: MB22_photospheric_ysol = 0.2485d0
290 : real(dp), parameter :: AAG21_photospheric_ysol = 0.2485d0
291 :
292 : character(len=iso_name_length) :: chem_element_main_iso_name(num_chem_elements)
293 : integer, parameter :: chem_element_name_len = iso_name_length
294 : character (len=chem_element_name_len) :: chem_element_Name(num_chem_elements)
295 : ! names for elements
296 :
297 :
298 : ! identifiers for different Z fractions.
299 : integer, parameter :: Custom_zfracs = 0
300 : integer, parameter :: AG89_zfracs = 1
301 : integer, parameter :: GN93_zfracs = 2
302 : integer, parameter :: GS98_zfracs = 3
303 : integer, parameter :: L03_zfracs = 4
304 : integer, parameter :: AGS05_zfracs = 5
305 : integer, parameter :: AGSS09_zfracs = 6
306 : integer, parameter :: L09_zfracs = 7
307 : integer, parameter :: A09_Prz_zfracs = 8
308 : integer, parameter :: MB22_photospheric_zfracs = 9
309 : integer, parameter :: AAG21_photospheric_zfracs = 10
310 :
311 :
312 : real(dp) :: AG89_element_zfrac(num_chem_elements) ! fraction by mass of total Z
313 : ! Anders & Grevesse 1989
314 :
315 : real(dp) :: GN93_element_zfrac(num_chem_elements) ! fraction by mass of total Z
316 : ! Grevesse and Noels 1993 abundances
317 :
318 : real(dp) :: GS98_element_zfrac(num_chem_elements) ! fraction by mass of total Z
319 : ! Grevesse and Sauval 1998 abundances
320 :
321 : real(dp) :: L03_element_zfrac(num_chem_elements) ! fraction by mass of total Z
322 : ! Lodders 2003 abundances
323 :
324 : real(dp) :: AGS05_element_zfrac(num_chem_elements) ! fraction by mass of total Z
325 : ! Asplund, Grevesse, and Sauval 2005 abundances
326 :
327 : real(dp) :: AGSS09_element_zfrac(num_chem_elements) ! fraction by mass of total Z
328 : ! Asplund, Grevesse, Sauval, and Scott 2009 abundances
329 : ! Annu. Rev. Astron. Astrophys. 2009. 47:481–522
330 :
331 : real(dp) :: L09_element_zfrac(num_chem_elements) ! fraction by mass of total Z
332 : ! Lodders and Palme, 2009. (http://adsabs.harvard.edu/abs/2009M%26PSA..72.5154L)
333 :
334 : real(dp) :: A09_Prz_zfrac(num_chem_elements) ! fraction by mass of the total Z
335 : ! Abundances are abased on Asplund, Grevesse, Sauval, and Scott 2009, ARA&A, 47:481–522
336 : ! but that of the some key elements are updated based on:
337 : ! "Present-day cosmic abundances ..."
338 : ! by Nieva, M.-F. & Przybilla, N. 2012, A&A, 539, 143
339 : ! and the proceeding paper
340 : ! "Hot Stars and Cosmic Abundances" by
341 : ! Przybilla N. Nieva M. F., Irrgang A. and Butler K. 2013, EAS Publ. Ser.
342 : ! in preparation
343 : ! The modified abundances w.r.t. A09 are (eps_El = log(El/H)+12.0)
344 : ! eps_He = 10.99
345 : ! eps_C = 8.33
346 : ! eps_N = 7.79
347 : ! eps_O = 8.76
348 : ! eps_Ne = 8.09 ! Cunha et al. (2006) give 8.11
349 : ! eps_Mg = 7.56
350 : ! eps_Al = 6.30
351 : ! eps_Si = 7.50
352 : ! eps_S = 7.14
353 : ! eps_Ar = 6.50
354 : ! eps_Fe = 7.52
355 :
356 : real(dp) :: MB22_photospheric_element_zfrac(num_chem_elements) ! fraction by mass of total Z
357 : ! Ekaterina Magg et al. , A&A 661, A140 (2022) photospheric abundance.
358 :
359 : real(dp) :: AAG21_photospheric_element_zfrac(num_chem_elements) ! fraction by mass of total Z
360 : ! Asplund et al. A&A 653, A141 (2021) photospheric abundance.
361 :
362 : type (integer_dict), pointer :: chem_element_names_dict
363 :
364 :
365 : real(dp) :: element_atomic_weight(num_chem_elements)
366 : ! de Laeter et al, Pure and Applied Chemistry 75(6), 683–799, 2003.
367 : ! (IUPAC Technical Report)
368 :
369 :
370 : ! temperature values at which partition function is defined
371 : real(dp), dimension(npart) :: Tpart
372 :
373 :
374 : ! mass excess of proton, neutron in MeV (for calculating binding energies)
375 : ! should be consistent with the mass excess of the prot and neut from the isotopes.data file
376 : ! Set in chem_isos_io.f90
377 : real(dp) :: del_Mp, del_Mn
378 :
379 : type nuclide_data
380 : integer :: nnuclides
381 : character(len=iso_name_length), dimension(:), pointer :: name ! name of nuclide
382 : integer, dimension(:), pointer :: chem_id ! (nnuclides)
383 : ! gives chem_id for member of nuclide_data
384 : integer, dimension(:), pointer :: nuclide ! (num_chem_isos)
385 : ! gives index in nuclide_data 1 .. nnuclides or 0 if not included
386 : real(dp), dimension(:), pointer :: W ! atomic weight (mass in amu units)
387 : integer, dimension(:), pointer :: Z ! number of protons
388 : integer, dimension(:), pointer :: N ! number of neutrons
389 : integer, dimension(:), pointer :: Z_plus_N ! number of baryons
390 : integer, dimension(:), pointer :: isomeric_state ! 0 is default
391 : real(dp), dimension(:), pointer :: spin ! ground-state spin
392 : real(dp), dimension(:), pointer :: binding_energy
393 : ! the binding energy is B = Z*del_Mp + N*del_Mn - mass_excess
394 : real(dp), dimension(:,:), pointer :: pfcn ! table of partition function
395 : real(dp), dimension(:), pointer :: mass_excess
396 : real(dp), dimension(:), pointer :: Z53 ! cache expensive Z^5/3 result
397 : end type nuclide_data
398 :
399 : type (nuclide_data) :: chem_isos ! from winvn
400 : type (integer_dict), pointer :: chem_isos_dict
401 : integer :: num_chem_isos ! no. entries in isotopes database
402 :
403 : ! storage container for a set of nuclides,
404 : ! used for extracting a subset of the full winvn database
405 : ! a set of nuclides is actually an array of these items
406 : type nuclide_set
407 : character(len=iso_name_length) :: nuclide
408 : integer :: rank
409 : end type nuclide_set
410 :
411 :
412 : ! reaction categories
413 :
414 : integer, parameter :: ipp = 1 ! pp chains
415 : integer, parameter :: icno = 2 ! cno cycles
416 : integer, parameter :: i3alf = 3 ! triple alpha
417 :
418 : ! "burn" in the following means decays or captures of protons, alphas, or neutrons
419 : integer, parameter :: i_burn_c = 4
420 : integer, parameter :: i_burn_n = 5
421 : integer, parameter :: i_burn_o = 6
422 : integer, parameter :: i_burn_ne = 7
423 : integer, parameter :: i_burn_na = 8
424 : integer, parameter :: i_burn_mg = 9
425 : integer, parameter :: i_burn_si = 10
426 : integer, parameter :: i_burn_s = 11
427 : integer, parameter :: i_burn_ar = 12
428 : integer, parameter :: i_burn_ca = 13
429 : integer, parameter :: i_burn_ti = 14
430 : integer, parameter :: i_burn_cr = 15
431 : integer, parameter :: i_burn_fe = 16
432 :
433 : integer, parameter :: icc = 17 ! c12 + c12
434 : integer, parameter :: ico = 18 ! c12 + o16
435 : integer, parameter :: ioo = 19 ! o16 + o16
436 :
437 : integer, parameter :: ipnhe4 = 20 ! 2prot + 2neut -> he4
438 :
439 : integer, parameter :: iphoto = 21 ! photodisintegration
440 : ! note: for photodisintegrations, eps_nuc will be negative.
441 :
442 : integer, parameter :: i_ni56_co56 = 22 ! ni56 -> co56
443 : integer, parameter :: i_co56_fe56 = 23 ! co56 -> fe56
444 : integer, parameter :: iother = 24 ! misc.
445 :
446 : integer, parameter :: num_categories = iother
447 :
448 : integer, parameter :: maxlen_category_name = 16
449 : character (len=maxlen_category_name) :: category_name(num_categories)
450 : type (integer_dict), pointer :: category_names_dict
451 :
452 :
453 : ! some commonly used values of get_nuclide_index
454 : integer :: &
455 : ih1, ih2, ih3, &
456 : ihe3, ihe4, &
457 : ili6, ili7, ili8, &
458 : ibe7, ibe8, ibe9, ibe10, ibe11, &
459 : ib8, ib10, ib11, ib12, ib13, ib14, &
460 : ic9, ic10, ic11, ic12, ic13, ic14, ic15, ic16, &
461 : in12, in13, in14, in15, in16, in17, in18, in19, in20, &
462 : io13, io14, io15, io16, io17, io18, io19, io20, &
463 : if15, if16, if17, if18, if19, if20, if21, if22, if23, if24, &
464 : ine17, ine18, ine19, ine20, ine21, ine22, ine23, ine24, ine25, &
465 : ine26, ine27, ine28, &
466 : ina20, ina21, ina22, ina23, ina24, ina25, ina26, ina27, ina28, &
467 : ina29, ina30, ina31, &
468 : img20, img21, img22, img23, img24, img25, img26, img27, img28, &
469 : img29, img30, img31, img32, img33, &
470 : ial22, ial23, ial24, ial25, ial26, ial27, ial28, ial29, ial30, &
471 : ial31, ial32, ial33, ial34, ial35, &
472 : isi22, isi23, isi24, isi25, isi26, isi27, isi28, isi29, isi30, &
473 : isi31, isi32, isi33, isi34, isi35, isi36, isi37, isi38, &
474 : ip26, ip27, ip28, ip29, ip30, ip31, ip32, ip33, ip34, ip35, &
475 : ip36, ip37, ip38, ip39, ip40, &
476 : is27, is28, is29, is30, is31, is32, is33, is34, is35, is36, &
477 : is37, is38, is39, is40, is41, is42, &
478 : icl31, icl32, icl33, icl34, icl35, icl36, icl37, icl38, icl39, &
479 : icl40, icl41, icl42, icl43, icl44, &
480 : iar31, iar32, iar33, iar34, iar35, iar36, iar37, iar38, iar39, &
481 : iar40, iar41, iar42, iar43, iar44, iar45, iar46, iar47, &
482 : ik35, ik36, ik37, ik38, ik39, ik40, ik41, ik42, ik43, ik44, ik45, ik46, ik47, &
483 : ica35, ica36, ica37, ica38, ica39, ica40, ica41, ica42, ica43, &
484 : ica44, ica45, ica46, ica47, ica48, ica49, ica50, ica51, ica52, ica53, &
485 : isc40, isc41, isc42, isc43, isc44, isc45, isc46, isc47, isc48, &
486 : isc49, isc50, isc51, isc52, isc53, &
487 : iti39, iti40, iti41, iti42, iti43, iti44, iti45, iti46, iti47, &
488 : iti48, iti49, iti50, iti51, iti52, iti53, iti54, iti55, &
489 : iv43, iv44, iv45, iv46, iv47, iv48, iv49, iv50, iv51, iv52, iv53, &
490 : iv54, iv55, iv56, iv57, &
491 : icr43, icr44, icr45, icr46, icr47, icr48, icr49, icr50, icr51, &
492 : icr52, icr53, icr54, icr55, icr56, icr57, icr58, icr59, &
493 : icr60, icr61, icr62, icr63, icr64, icr65, icr66, &
494 : imn46, imn47, imn48, imn49, imn50, imn51, imn52, imn53, imn54, &
495 : imn55, imn56, imn57, imn58, imn59, imn60, imn61, imn62, imn63, &
496 : ife46, ife47, ife48, ife49, ife50, ife51, ife52, ife53, ife54, ife55, ife56, &
497 : ife57, ife58, ife59, ife60, ife61, ife62, ife63, ife64, ife65, ife66, ife68, &
498 : ico50, ico51, ico52, ico53, ico54, ico55, ico56, ico57, ico58, &
499 : ico59, ico60, ico61, ico62, ico63, ico64, ico65, ico66, ico67, &
500 : ini50, ini51, ini52, ini53, ini54, ini55, ini56, ini57, ini58, &
501 : ini59, ini60, ini61, ini62, ini63, ini64, ini65, ini66, ini67, &
502 : ini68, ini69, ini70, ini71, ini72, ini73, &
503 : icu56, icu57, icu58, icu59, icu60, icu61, icu62, icu63, icu64, &
504 : icu65, icu66, icu67, icu68, icu69, icu70, icu71, icu72, &
505 : izn55, izn56, izn57, izn58, izn59, izn60, izn61, izn62, izn63, izn64, &
506 : izn65, izn66, izn67, izn68, izn69, izn70, izn71, izn72, izn73, izn74, &
507 : iga60, iga61, iga62, iga63, iga64, iga65, iga66, iga67, iga68, &
508 : iga69, iga70, iga71, iga72, iga73, iga74, iga75, &
509 : ige59, ige60, ige61, ige62, ige63, ige64, ige65, ige66, ige67, &
510 : ige68, ige69, ige70, ige71, ige72, ige73, ige74, ige75, ige76, &
511 : ias71,ias72, ias73, ias74, ias75, ias76, ias77, ias78, ias79, &
512 : ise68, ise69, ise70, ise71, ise72, ise73, ise74, ise75, ise76, &
513 : ikr70, ikr72, isr74, isr75, isr76, izr77, izr80, imo82, imo84, &
514 : iru86, iru87, iru88, ipd89, ipd91, ipd92, icd93, icd96, &
515 : isn98, isn100, isn102, isn104, &
516 : ineut, iprot
517 :
518 :
519 : logical :: chem_has_been_initialized = .false.
520 :
521 :
522 : contains
523 :
524 :
525 16 : subroutine init_chem_tables
526 : use utils_lib, only: integer_dict_define, integer_dict_create_hash
527 :
528 : integer :: i, ierr
529 :
530 : Tpart = [ &
531 : 0.10d0, 0.15d0, 0.20d0, 0.30d0, 0.40d0, 0.50d0, &
532 : 0.60d0, 0.70d0, 0.80d0, 0.90d0, 1.00d0, 1.50d0, &
533 : 2.00d0, 2.50d0, 3.00d0, 3.50d0, 4.00d0, 4.50d0, &
534 16 : 5.00d0, 6.00d0, 7.00d0, 8.00d0, 9.00d0, 10.0d0 ]
535 :
536 16 : call init_ag_data
537 :
538 16 : call init_chem_element_names
539 :
540 16 : call init_chem_element_main_iso_names
541 :
542 16 : call init_element_atomic_weights
543 :
544 16 : call init_AG89_data
545 :
546 16 : call init_GN93_data
547 :
548 16 : call init_GS98_data
549 :
550 16 : call init_L03_data
551 :
552 16 : call init_AGS05_data
553 :
554 16 : call init_AGSS09_data
555 :
556 16 : call init_A09_Przybilla_data
557 :
558 16 : call init_MB22_photospheric_data
559 :
560 16 : call init_AAG21_photospheric_data
561 :
562 16 : call init_L09_data
563 :
564 16 : nullify(chem_element_names_dict)
565 1808 : do i=1,num_chem_elements
566 1792 : call integer_dict_define(chem_element_names_dict, chem_element_Name(i), i, ierr)
567 1808 : if (ierr /= 0) then
568 0 : write(*,*) 'FATAL ERROR: init_chem_tables failed in integer_dict_define'
569 0 : flush(6)
570 0 : call mesa_error(__FILE__,__LINE__)
571 : end if
572 : end do
573 16 : call integer_dict_create_hash(chem_element_names_dict, ierr)
574 16 : if (ierr /= 0) then
575 0 : write(*,*) 'FATAL ERROR: init_chem_tables failed in integer_dict_create_hash'
576 0 : flush(6)
577 0 : call mesa_error(__FILE__,__LINE__)
578 : end if
579 :
580 16 : call set_category_names
581 16 : nullify(category_names_dict)
582 400 : do i=1,num_categories
583 384 : call integer_dict_define(category_names_dict, category_name(i), i, ierr)
584 400 : if (ierr /= 0) then
585 0 : write(*,*) 'FATAL ERROR: rates_def_init failed in integer_dict_define'
586 0 : flush(6)
587 0 : call mesa_error(__FILE__,__LINE__)
588 : end if
589 : end do
590 16 : call integer_dict_create_hash(category_names_dict, ierr)
591 16 : if (ierr /= 0) then
592 0 : write(*,*) 'FATAL ERROR: rates_def_init failed in integer_dict_create_hash'
593 0 : flush(6)
594 0 : call mesa_error(__FILE__,__LINE__)
595 : end if
596 :
597 32 : end subroutine init_chem_tables
598 :
599 :
600 16 : subroutine init_ag_data
601 : use utils_lib
602 :
603 16 : real(dp) :: sum
604 : integer :: i, j, ierr
605 :
606 : !..names of the stable isotopes
607 : namsol(1:120) = [ &
608 : 'h1 ','h2 ','he3 ','he4 ','li6 ','li7 ','be9 ','b10 ', &
609 : 'b11 ','c12 ','c13 ','n14 ','n15 ','o16 ','o17 ','o18 ', &
610 : 'f19 ','ne20 ','ne21 ','ne22 ','na23 ','mg24 ','mg25 ','mg26 ', &
611 : 'al27 ','si28 ','si29 ','si30 ','p31 ','s32 ','s33 ','s34 ', &
612 : 's36 ','cl35 ','cl37 ','ar36 ','ar38 ','ar40 ','k39 ','k40 ', &
613 : 'k41 ','ca40 ','ca42 ','ca43 ','ca44 ','ca46 ','ca48 ','sc45 ', &
614 : 'ti46 ','ti47 ','ti48 ','ti49 ','ti50 ','v50 ','v51 ','cr50 ', &
615 : 'cr52 ','cr53 ','cr54 ','mn55 ','fe54 ','fe56 ','fe57 ','fe58 ', &
616 : 'co59 ','ni58 ','ni60 ','ni61 ','ni62 ','ni64 ','cu63 ','cu65 ', &
617 : 'zn64 ','zn66 ','zn67 ','zn68 ','zn70 ','ga69 ','ga71 ','ge70 ', &
618 : 'ge72 ','ge73 ','ge74 ','ge76 ','as75 ','se74 ','se76 ','se77 ', &
619 : 'se78 ','se80 ','se82 ','br79 ','br81 ','kr78 ','kr80 ','kr82 ', &
620 : 'kr83 ','kr84 ','kr86 ','rb85 ','rb87 ','sr84 ','sr86 ','sr87 ', &
621 : 'sr88 ','y89 ','zr90 ','zr91 ','zr92 ','zr94 ','zr96 ','nb93 ', &
622 1936 : 'mo92 ','mo94 ','mo95 ','mo96 ','mo97 ','mo98 ','mo100','ru96 ' ]
623 :
624 : namsol(121:240) = [ &
625 : 'ru98 ','ru99 ','ru100','ru101','ru102','ru104','rh103','pd102', &
626 : 'pd104','pd105','pd106','pd108','pd110','ag107','ag109','cd106', &
627 : 'cd108','cd110','cd111','cd112','cd113','cd114','cd116','in113', &
628 : 'in115','sn112','sn114','sn115','sn116','sn117','sn118','sn119', &
629 : 'sn120','sn122','sn124','sb121','sb123','te120','te122','te123', &
630 : 'te124','te125','te126','te128','te130','i127 ','xe124','xe126', &
631 : 'xe128','xe129','xe130','xe131','xe132','xe134','xe136','cs133', &
632 : 'ba130','ba132','ba134','ba135','ba136','ba137','ba138','la138', &
633 : 'la139','ce136','ce138','ce140','ce142','pr141','nd142','nd143', &
634 : 'nd144','nd145','nd146','nd148','nd150','sm144','sm147','sm148', &
635 : 'sm149','sm150','sm152','sm154','eu151','eu153','gd152','gd154', &
636 : 'gd155','gd156','gd157','gd158','gd160','tb159','dy156','dy158', &
637 : 'dy160','dy161','dy162','dy163','dy164','ho165','er162','er164', &
638 : 'er166','er167','er168','er170','tm169','yb168','yb170','yb171', &
639 1936 : 'yb172','yb173','yb174','yb176','lu175','lu176','hf174','hf176' ]
640 :
641 : namsol(241:286) = [ &
642 : 'hf177','hf178','hf179','hf180','ta180','ta181','w180 ','w182 ', &
643 : 'w183 ','w184 ','w186 ','re185','re187','os184','os186','os187', &
644 : 'os188','os189','os190','os192','ir191','ir193','pt190','pt192', &
645 : 'pt194','pt195','pt196','pt198','au197','hg196','hg198','hg199', &
646 : 'hg200','hg201','hg202','hg204','tl203','tl205','pb204','pb206', &
647 752 : 'pb207','pb208','bi209','th232','u235 ','u238 ' ]
648 :
649 :
650 : !..anders & grevesse 1989 solar mass fractions
651 : solx(1:45) = [ &
652 : 7.0573D-01, 4.8010D-05, 2.9291D-05, 2.7521D-01, 6.4957D-10, &
653 : 9.3490D-09, 1.6619D-10, 1.0674D-09, 4.7301D-09, 3.0324D-03, &
654 : 3.6501D-05, 1.1049D-03, 4.3634D-06, 9.5918D-03, 3.8873D-06, &
655 : 2.1673D-05, 4.0515D-07, 1.6189D-03, 4.1274D-06, 1.3022D-04, &
656 : 3.3394D-05, 5.1480D-04, 6.7664D-05, 7.7605D-05, 5.8052D-05, &
657 : 6.5301D-04, 3.4257D-05, 2.3524D-05, 8.1551D-06, 3.9581D-04, &
658 : 3.2221D-06, 1.8663D-05, 9.3793D-08, 2.5320D-06, 8.5449D-07, &
659 : 7.7402D-05, 1.5379D-05, 2.6307D-08, 3.4725D-06, 4.4519D-10, &
660 736 : 2.6342D-07, 5.9898D-05, 4.1964D-07, 8.9734D-07, 1.4135D-06 ]
661 :
662 : solx(46:90) = [ &
663 : 2.7926D-09, 1.3841D-07, 3.8929D-08, 2.2340D-07, 2.0805D-07, &
664 : 2.1491D-06, 1.6361D-07, 1.6442D-07, 9.2579D-10, 3.7669D-07, &
665 : 7.4240D-07, 1.4863D-05, 1.7160D-06, 4.3573D-07, 1.3286D-05, &
666 : 7.1301D-05, 1.1686D-03, 2.8548D-05, 3.6971D-06, 3.3579D-06, &
667 : 4.9441D-05, 1.9578D-05, 8.5944D-07, 2.7759D-06, 7.2687D-07, &
668 : 5.7528D-07, 2.6471D-07, 9.9237D-07, 5.8765D-07, 8.7619D-08, &
669 : 4.0593D-07, 1.3811D-08, 3.9619D-08, 2.7119D-08, 4.3204D-08, &
670 : 5.9372D-08, 1.7136D-08, 8.1237D-08, 1.7840D-08, 1.2445D-08, &
671 736 : 1.0295D-09, 1.0766D-08, 9.1542D-09, 2.9003D-08, 6.2529D-08 ]
672 :
673 : solx(91:135) = [ &
674 : 1.1823D-08, 1.1950D-08, 1.2006D-08, 3.0187D-10, 2.0216D-09, &
675 : 1.0682D-08, 1.0833D-08, 5.4607D-08, 1.7055D-08, 1.1008D-08, &
676 : 4.3353D-09, 2.8047D-10, 5.0468D-09, 3.6091D-09, 4.3183D-08, &
677 : 1.0446D-08, 1.3363D-08, 2.9463D-09, 4.5612D-09, 4.7079D-09, &
678 : 7.7706D-10, 1.6420D-09, 8.7966D-10, 5.6114D-10, 9.7562D-10, &
679 : 1.0320D-09, 5.9868D-10, 1.5245D-09, 6.2225D-10, 2.5012D-10, &
680 : 8.6761D-11, 5.9099D-10, 5.9190D-10, 8.0731D-10, 1.5171D-09, &
681 : 9.1547D-10, 8.9625D-10, 3.6637D-11, 4.0775D-10, 8.2335D-10, &
682 736 : 1.0189D-09, 1.0053D-09, 4.5354D-10, 6.8205D-10, 6.4517D-10 ]
683 :
684 : solx(136:180) = [ &
685 : 5.3893D-11, 3.9065D-11, 5.5927D-10, 5.7839D-10, 1.0992D-09, &
686 : 5.6309D-10, 1.3351D-09, 3.5504D-10, 2.2581D-11, 5.1197D-10, &
687 : 1.0539D-10, 7.1802D-11, 3.9852D-11, 1.6285D-09, 8.6713D-10, &
688 : 2.7609D-09, 9.8731D-10, 3.7639D-09, 5.4622D-10, 6.9318D-10, &
689 : 5.4174D-10, 4.1069D-10, 1.3052D-11, 3.8266D-10, 1.3316D-10, &
690 : 7.1827D-10, 1.0814D-09, 3.1553D-09, 4.9538D-09, 5.3600D-09, &
691 : 2.8912D-09, 1.7910D-11, 1.6223D-11, 3.3349D-10, 4.1767D-09, &
692 : 6.7411D-10, 3.3799D-09, 4.1403D-09, 1.5558D-09, 1.2832D-09, &
693 736 : 1.2515D-09, 1.5652D-11, 1.5125D-11, 3.6946D-10, 1.0108D-09 ]
694 :
695 : solx(181:225) = [ &
696 : 1.2144D-09, 1.7466D-09, 1.1240D-08, 1.3858D-12, 1.5681D-09, &
697 : 7.4306D-12, 9.9136D-12, 3.5767D-09, 4.5258D-10, 5.9562D-10, &
698 : 8.0817D-10, 3.6533D-10, 7.1757D-10, 2.5198D-10, 5.2441D-10, &
699 : 1.7857D-10, 1.7719D-10, 2.9140D-11, 1.4390D-10, 1.0931D-10, &
700 : 1.3417D-10, 7.2470D-11, 2.6491D-10, 2.2827D-10, 1.7761D-10, &
701 : 1.9660D-10, 2.5376D-12, 2.8008D-11, 1.9133D-10, 2.6675D-10, &
702 : 2.0492D-10, 3.2772D-10, 2.9180D-10, 2.8274D-10, 8.6812D-13, &
703 : 1.4787D-12, 3.7315D-11, 3.0340D-10, 4.1387D-10, 4.0489D-10, &
704 736 : 4.6047D-10, 3.7104D-10, 1.4342D-12, 1.6759D-11, 3.5397D-10 ]
705 :
706 : solx(226:270) = [ &
707 : 2.4332D-10, 2.8557D-10, 1.6082D-10, 1.6159D-10, 1.3599D-12, &
708 : 3.2509D-11, 1.5312D-10, 2.3624D-10, 1.7504D-10, 3.4682D-10, &
709 : 1.4023D-10, 1.5803D-10, 4.2293D-12, 1.0783D-12, 3.4992D-11, &
710 : 1.2581D-10, 1.8550D-10, 9.3272D-11, 2.4131D-10, 1.1292D-14, &
711 : 9.4772D-11, 7.8768D-13, 1.6113D-10, 8.7950D-11, 1.8989D-10, &
712 : 1.7878D-10, 9.0315D-11, 1.5326D-10, 5.6782D-13, 5.0342D-11, &
713 : 5.1086D-11, 4.2704D-10, 5.2110D-10, 8.5547D-10, 1.3453D-09, &
714 : 1.1933D-09, 2.0211D-09, 8.1702D-13, 5.0994D-11, 2.1641D-09, &
715 736 : 2.2344D-09, 1.6757D-09, 4.8231D-10, 9.3184D-10, 2.3797D-12 ]
716 :
717 : solx(271:286) = [ &
718 : 1.7079D-10, 2.8843D-10, 3.9764D-10, 2.2828D-10, 5.1607D-10, &
719 : 1.2023D-10, 2.7882D-10, 6.7411D-10, 3.1529D-10, 3.1369D-09, &
720 : 3.4034D-09, 9.6809D-09, 7.6127D-10, 1.9659D-10, 3.8519D-13, &
721 272 : 5.3760D-11 ]
722 :
723 :
724 : !..charge of the stable isotopes
725 :
726 : izsol(1:117) = [ &
727 : 1, 1, 2, 2, 3, 3, 4, 5, 5, 6, 6, 7, 7, &
728 : 8, 8, 8, 9, 10, 10, 10, 11, 12, 12, 12, 13, 14, &
729 : 14, 14, 15, 16, 16, 16, 16, 17, 17, 18, 18, 18, 19, &
730 : 19, 19, 20, 20, 20, 20, 20, 20, 21, 22, 22, 22, 22, &
731 : 22, 23, 23, 24, 24, 24, 24, 25, 26, 26, 26, 26, 27, &
732 : 28, 28, 28, 28, 28, 29, 29, 30, 30, 30, 30, 30, 31, &
733 : 31, 32, 32, 32, 32, 32, 33, 34, 34, 34, 34, 34, 34, &
734 : 35, 35, 36, 36, 36, 36, 36, 36, 37, 37, 38, 38, 38, &
735 1888 : 38, 39, 40, 40, 40, 40, 40, 41, 42, 42, 42, 42, 42 ]
736 :
737 : izsol(118:234) = [ &
738 : 42, 42, 44, 44, 44, 44, 44, 44, 44, 45, 46, 46, 46, &
739 : 46, 46, 46, 47, 47, 48, 48, 48, 48, 48, 48, 48, 48, &
740 : 49, 49, 50, 50, 50, 50, 50, 50, 50, 50, 50, 50, 51, &
741 : 51, 52, 52, 52, 52, 52, 52, 52, 52, 53, 54, 54, 54, &
742 : 54, 54, 54, 54, 54, 54, 55, 56, 56, 56, 56, 56, 56, &
743 : 56, 57, 57, 58, 58, 58, 58, 59, 60, 60, 60, 60, 60, &
744 : 60, 60, 62, 62, 62, 62, 62, 62, 62, 63, 63, 64, 64, &
745 : 64, 64, 64, 64, 64, 65, 66, 66, 66, 66, 66, 66, 66, &
746 1888 : 67, 68, 68, 68, 68, 68, 68, 69, 70, 70, 70, 70, 70 ]
747 :
748 : izsol(235:286) = [ &
749 : 70, 70, 71, 71, 72, 72, 72, 72, 72, 72, 73, 73, 74, &
750 : 74, 74, 74, 74, 75, 75, 76, 76, 76, 76, 76, 76, 76, &
751 : 77, 77, 78, 78, 78, 78, 78, 78, 79, 80, 80, 80, 80, &
752 848 : 80, 80, 80, 81, 81, 82, 82, 82, 82, 83, 90, 92, 92 ]
753 :
754 :
755 : !..number of nucleons (protons and neutrons) in the stable isotopes
756 :
757 : iasol(1:117) = [ &
758 : 1, 2, 3, 4, 6, 7, 9, 10, 11, 12, 13, 14, 15, &
759 : 16, 17, 18, 19, 20, 21, 22, 23, 24, 25, 26, 27, 28, &
760 : 29, 30, 31, 32, 33, 34, 36, 35, 37, 36, 38, 40, 39, &
761 : 40, 41, 40, 42, 43, 44, 46, 48, 45, 46, 47, 48, 49, &
762 : 50, 50, 51, 50, 52, 53, 54, 55, 54, 56, 57, 58, 59, &
763 : 58, 60, 61, 62, 64, 63, 65, 64, 66, 67, 68, 70, 69, &
764 : 71, 70, 72, 73, 74, 76, 75, 74, 76, 77, 78, 80, 82, &
765 : 79, 81, 78, 80, 82, 83, 84, 86, 85, 87, 84, 86, 87, &
766 1888 : 88, 89, 90, 91, 92, 94, 96, 93, 92, 94, 95, 96, 97 ]
767 :
768 : iasol(118:234) = [ &
769 : 98, 100, 96, 98, 99, 100, 101, 102, 104, 103, 102, 104, 105, &
770 : 106, 108, 110, 107, 109, 106, 108, 110, 111, 112, 113, 114, 116, &
771 : 113, 115, 112, 114, 115, 116, 117, 118, 119, 120, 122, 124, 121, &
772 : 123, 120, 122, 123, 124, 125, 126, 128, 130, 127, 124, 126, 128, &
773 : 129, 130, 131, 132, 134, 136, 133, 130, 132, 134, 135, 136, 137, &
774 : 138, 138, 139, 136, 138, 140, 142, 141, 142, 143, 144, 145, 146, &
775 : 148, 150, 144, 147, 148, 149, 150, 152, 154, 151, 153, 152, 154, &
776 : 155, 156, 157, 158, 160, 159, 156, 158, 160, 161, 162, 163, 164, &
777 1888 : 165, 162, 164, 166, 167, 168, 170, 169, 168, 170, 171, 172, 173 ]
778 :
779 : iasol(235:286) = [ &
780 : 174, 176, 175, 176, 174, 176, 177, 178, 179, 180, 180, 181, 180, &
781 : 182, 183, 184, 186, 185, 187, 184, 186, 187, 188, 189, 190, 192, &
782 : 191, 193, 190, 192, 194, 195, 196, 198, 197, 196, 198, 199, 200, &
783 848 : 201, 202, 204, 203, 205, 204, 206, 207, 208, 209, 232, 235, 238 ]
784 :
785 :
786 : ! jcode tells the type progenitors each stable species can have.
787 : ! jcode = 0 if the species is the only stable one of that a
788 : ! = 1 if the species can have proton-rich progenitors
789 : ! = 2 if the species can have neutron-rich progenitors
790 : ! = 3 if the species can only be made as itself (eg k40)
791 :
792 : jcode(1:286) = [ &
793 : 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
794 : 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, &
795 : 0, 0, 0, 0, 0, 0, 2, 0, 0, 1, 0, 2, 0, &
796 : 3, 0, 1, 0, 0, 0, 2, 2, 0, 1, 0, 1, 0, &
797 : 2, 3, 0, 1, 0, 0, 2, 0, 1, 0, 0, 2, 0, &
798 : 1, 0, 0, 0, 2, 0, 0, 1, 0, 0, 0, 2, 0, &
799 : 0, 1, 0, 0, 2, 2, 0, 1, 1, 0, 2, 2, 2, &
800 : 0, 0, 1, 1, 1, 0, 2, 2, 0, 2, 1, 1, 1, &
801 : 0, 0, 0, 0, 2, 2, 2, 0, 1, 1, 0, 3, 0, &
802 : 2, 2, 1, 1, 0, 1, 0, 2, 2, 0, 1, 1, 0, &
803 : 2, 2, 2, 0, 0, 1, 1, 1, 0, 2, 2, 2, 2, &
804 : 1, 2, 1, 1, 1, 1, 0, 0, 0, 2, 2, 2, 0, &
805 : 2, 1, 1, 1, 3, 0, 2, 2, 2, 0, 1, 1, 1, &
806 : 0, 3, 0, 2, 2, 2, 0, 1, 1, 1, 0, 3, 0, &
807 : 2, 3, 0, 1, 1, 0, 2, 0, 1, 0, 2, 0, 0, &
808 : 2, 2, 1, 0, 1, 0, 1, 2, 2, 0, 0, 1, 1, &
809 : 0, 2, 0, 2, 2, 0, 1, 1, 1, 0, 2, 0, 2, &
810 : 0, 1, 1, 0, 0, 2, 2, 0, 1, 1, 0, 0, 0, &
811 : 2, 2, 0, 3, 1, 1, 0, 0, 0, 2, 3, 0, 1, &
812 : 0, 0, 2, 2, 0, 2, 1, 1, 1, 0, 0, 2, 2, &
813 : 0, 0, 1, 1, 0, 0, 2, 2, 0, 1, 1, 0, 0, &
814 16 : 0, 0, 2, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0]
815 :
816 :
817 : ! get sum and stuff residual into hydrogen
818 16 : sum = 0.0d0
819 4592 : do j=1,solsiz
820 4592 : sum = sum + solx(j)
821 : end do
822 16 : sum = 1.0d0 - sum
823 16 : solx(1) = solx(1) + sum
824 :
825 16 : sum = 0.0d0
826 4592 : do j=1,solsiz
827 4592 : if (izsol(j) >= 3) then
828 4512 : sum = sum + solx(j)
829 : end if
830 : end do
831 16 : zsol = sum
832 :
833 16 : sum = 0.0d0
834 4592 : do j=1,solsiz
835 4592 : sum = sum + dble(izsol(j))/dble(iasol(j))*solx(j)
836 : end do
837 16 : yesol = sum
838 :
839 16 : nullify(Xsol_names_dict)
840 : ierr = 0
841 4592 : do i=1,solsiz
842 4576 : call integer_dict_define(Xsol_names_dict, namsol(i), i, ierr)
843 4592 : if (ierr /= 0) then
844 0 : write(*,*) 'FATAL ERROR: chem_init, init_ag_data failed in integer_dict_define'
845 0 : flush(6)
846 0 : call mesa_error(__FILE__,__LINE__)
847 : end if
848 : end do
849 16 : call integer_dict_define(Xsol_names_dict, 'p', 1, ierr)
850 16 : if (ierr /= 0) then
851 0 : write(*,*) 'FATAL ERROR: chem_init, init_ag_data failed in integer_dict_define'
852 0 : flush(6)
853 0 : call mesa_error(__FILE__,__LINE__)
854 : end if
855 16 : call integer_dict_define(Xsol_names_dict, 'd', 2, ierr)
856 16 : if (ierr /= 0) then
857 0 : write(*,*) 'FATAL ERROR: chem_init, init_ag_data failed in integer_dict_define'
858 0 : flush(6)
859 0 : call mesa_error(__FILE__,__LINE__)
860 : end if
861 16 : call integer_dict_create_hash(Xsol_names_dict, ierr)
862 16 : if (ierr /= 0) then
863 0 : write(*,*) 'FATAL ERROR: chem_init, init_ag_data failed in integer_dict_create_hash'
864 0 : flush(6)
865 0 : call mesa_error(__FILE__,__LINE__)
866 : end if
867 :
868 16 : end subroutine init_ag_data
869 :
870 :
871 16 : subroutine init_chem_element_names
872 : integer :: i
873 1808 : chem_element_name(:) = ''
874 :
875 1808 : chem_element_name(1:num_chem_elements) = el_name(1:num_chem_elements)
876 :
877 1808 : do i=1,num_chem_elements
878 1808 : if (len_trim(chem_element_name(i)) == 0) then
879 0 : write(*,*)'missing chem_element_name(i)', i
880 0 : flush(6)
881 0 : call mesa_error(__FILE__,__LINE__)
882 : end if
883 : end do
884 :
885 16 : end subroutine init_chem_element_names
886 :
887 :
888 16 : subroutine init_chem_element_main_iso_names
889 : ! the iso with the largest number abundance according to Lodders03
890 : integer :: i
891 1808 : chem_element_main_iso_name(:) = ''
892 :
893 : !periodic table, row 1
894 16 : chem_element_main_iso_name(e_h) = 'h1'
895 16 : chem_element_main_iso_name(e_he) = 'he4'
896 :
897 : !periodic table, row 2
898 16 : chem_element_main_iso_name(e_li) = 'li7'
899 16 : chem_element_main_iso_name(e_be) = 'be9'
900 16 : chem_element_main_iso_name(e_b) = 'b11'
901 16 : chem_element_main_iso_name(e_c) = 'c12'
902 16 : chem_element_main_iso_name(e_n) = 'n14'
903 16 : chem_element_main_iso_name(e_o) = 'o16'
904 16 : chem_element_main_iso_name(e_f) = 'f19'
905 16 : chem_element_main_iso_name(e_ne) = 'ne20'
906 :
907 : !periodic table, row 3
908 16 : chem_element_main_iso_name(e_na) = 'na23'
909 16 : chem_element_main_iso_name(e_mg) = 'mg24'
910 16 : chem_element_main_iso_name(e_al) = 'al27'
911 16 : chem_element_main_iso_name(e_si) = 'si28'
912 16 : chem_element_main_iso_name(e_p) = 'p31'
913 16 : chem_element_main_iso_name(e_s) = 's32'
914 16 : chem_element_main_iso_name(e_cl) = 'cl35'
915 16 : chem_element_main_iso_name(e_ar) = 'ar36'
916 :
917 : !periodic table, row 4
918 16 : chem_element_main_iso_name(e_k) = 'k39'
919 16 : chem_element_main_iso_name(e_ca) = 'ca40'
920 16 : chem_element_main_iso_name(e_sc) = 'sc45'
921 16 : chem_element_main_iso_name(e_ti) = 'ti48'
922 16 : chem_element_main_iso_name(e_v) = 'v51'
923 16 : chem_element_main_iso_name(e_cr) = 'cr52'
924 16 : chem_element_main_iso_name(e_mn) = 'mn55'
925 16 : chem_element_main_iso_name(e_fe) = 'fe56'
926 16 : chem_element_main_iso_name(e_co) = 'co59'
927 16 : chem_element_main_iso_name(e_ni) = 'ni58'
928 16 : chem_element_main_iso_name(e_cu) = 'cu63'
929 16 : chem_element_main_iso_name(e_zn) = 'zn64'
930 16 : chem_element_main_iso_name(e_ga) = 'ga69'
931 16 : chem_element_main_iso_name(e_ge) = 'ge74'
932 16 : chem_element_main_iso_name(e_as) = 'as75'
933 16 : chem_element_main_iso_name(e_se) = 'se80'
934 16 : chem_element_main_iso_name(e_br) = 'br79'
935 16 : chem_element_main_iso_name(e_kr) = 'kr84'
936 :
937 : !periodic table, row 5
938 16 : chem_element_main_iso_name(e_rb) = 'rb85'
939 16 : chem_element_main_iso_name(e_sr) = 'sr88'
940 16 : chem_element_main_iso_name(e_y) = 'y89'
941 16 : chem_element_main_iso_name(e_zr) = 'zr90'
942 16 : chem_element_main_iso_name(e_nb) = 'nb93'
943 16 : chem_element_main_iso_name(e_mo) = 'mo98'
944 16 : chem_element_main_iso_name(e_tc) = 'tc97'
945 16 : chem_element_main_iso_name(e_ru) = 'ru102'
946 16 : chem_element_main_iso_name(e_rh) = 'rh103'
947 16 : chem_element_main_iso_name(e_pd) = 'pd106'
948 16 : chem_element_main_iso_name(e_ag) = 'ag107'
949 16 : chem_element_main_iso_name(e_cd) = 'cd114'
950 16 : chem_element_main_iso_name(e_in) = 'in115'
951 16 : chem_element_main_iso_name(e_sn) = 'sn120'
952 16 : chem_element_main_iso_name(e_sb) = 'sb121'
953 16 : chem_element_main_iso_name(e_te) = 'te130'
954 16 : chem_element_main_iso_name(e_i ) = 'i127'
955 16 : chem_element_main_iso_name(e_xe) = 'xe132'
956 :
957 : !periodic table, row 6
958 16 : chem_element_main_iso_name(e_cs) = 'cs133'
959 16 : chem_element_main_iso_name(e_ba) = 'ba138'
960 16 : chem_element_main_iso_name(e_la) = 'la139'
961 16 : chem_element_main_iso_name(e_ce) = 'ce140'
962 16 : chem_element_main_iso_name(e_pr) = 'pr141'
963 16 : chem_element_main_iso_name(e_nd) = 'nd142'
964 16 : chem_element_main_iso_name(e_pm) = 'pm145'
965 16 : chem_element_main_iso_name(e_sm) = 'sm152'
966 16 : chem_element_main_iso_name(e_eu) = 'eu153'
967 16 : chem_element_main_iso_name(e_gd) = 'gd158'
968 16 : chem_element_main_iso_name(e_tb) = 'tb159'
969 16 : chem_element_main_iso_name(e_dy) = 'dy164'
970 16 : chem_element_main_iso_name(e_ho) = 'ho165'
971 16 : chem_element_main_iso_name(e_er) = 'er166'
972 16 : chem_element_main_iso_name(e_tm) = 'tm169'
973 16 : chem_element_main_iso_name(e_yb) = 'yb174'
974 16 : chem_element_main_iso_name(e_lu) = 'lu175'
975 16 : chem_element_main_iso_name(e_hf) = 'hf180'
976 16 : chem_element_main_iso_name(e_ta) = 'ta181'
977 16 : chem_element_main_iso_name(e_w ) = 'w184'
978 16 : chem_element_main_iso_name(e_re) = 're187'
979 16 : chem_element_main_iso_name(e_os) = 'os192'
980 16 : chem_element_main_iso_name(e_ir) = 'ir193'
981 16 : chem_element_main_iso_name(e_pt) = 'pt195'
982 16 : chem_element_main_iso_name(e_au) = 'au197'
983 16 : chem_element_main_iso_name(e_hg) = 'hg202'
984 16 : chem_element_main_iso_name(e_tl) = 'tl205'
985 16 : chem_element_main_iso_name(e_pb) = 'pb208'
986 16 : chem_element_main_iso_name(e_bi) = 'bi209'
987 16 : chem_element_main_iso_name(e_po) = 'po210'
988 16 : chem_element_main_iso_name(e_at) = 'at210'
989 : ! need a rule here to continue for the unstable nuclei -- longest half-life?
990 16 : chem_element_main_iso_name(e_rn) = 'rn222'
991 16 : chem_element_main_iso_name(e_fr) = 'fr223'
992 16 : chem_element_main_iso_name(e_ra) = 'ra226'
993 16 : chem_element_main_iso_name(e_ac) = 'ac227'
994 16 : chem_element_main_iso_name(e_th) = 'th232'
995 16 : chem_element_main_iso_name(e_pa) = 'pa231'
996 16 : chem_element_main_iso_name(e_u) = 'u238'
997 16 : chem_element_main_iso_name(e_np) = 'np237'
998 16 : chem_element_main_iso_name(e_pu) = 'pu244'
999 16 : chem_element_main_iso_name(e_am) = 'am243'
1000 16 : chem_element_main_iso_name(e_cm) = 'cm247'
1001 16 : chem_element_main_iso_name(e_bk) = 'bk247'
1002 16 : chem_element_main_iso_name(e_cf) = 'cf251'
1003 16 : chem_element_main_iso_name(e_es) = 'es252'
1004 16 : chem_element_main_iso_name(e_fm) = 'fm257'
1005 16 : chem_element_main_iso_name(e_md) = 'md258'
1006 16 : chem_element_main_iso_name(e_no) = 'no259'
1007 16 : chem_element_main_iso_name(e_lr) = 'lr262'
1008 16 : chem_element_main_iso_name(e_rf) = 'rf261'
1009 16 : chem_element_main_iso_name(e_db) = 'db268'
1010 16 : chem_element_main_iso_name(e_sg) = 'sg271'
1011 16 : chem_element_main_iso_name(e_bh) = 'bh274'
1012 16 : chem_element_main_iso_name(e_hs) = 'hs270'
1013 16 : chem_element_main_iso_name(e_mt) = 'mt278'
1014 16 : chem_element_main_iso_name(e_ds) = 'ds281'
1015 16 : chem_element_main_iso_name(e_rg) = 'rg281'
1016 16 : chem_element_main_iso_name(e_cn) = 'cn285'
1017 :
1018 :
1019 1808 : do i=1,num_chem_elements
1020 1808 : if (len_trim(chem_element_main_iso_name(i)) == 0) then
1021 0 : write(*,*)'missing chem_element_main_iso_name', i
1022 0 : flush(6)
1023 0 : call mesa_error(__FILE__,__LINE__)
1024 : end if
1025 : end do
1026 :
1027 16 : end subroutine init_chem_element_main_iso_names
1028 :
1029 :
1030 16 : subroutine init_element_atomic_weights
1031 : use utils_lib, only: integer_dict_lookup
1032 : integer :: i, isotope_index, ierr
1033 16 : element_atomic_weight(:) = 0d0
1034 :
1035 1808 : do i = e_h, e_cn
1036 1792 : call integer_dict_lookup(chem_isos_dict, chem_element_main_iso_name(i), isotope_index, ierr)
1037 1792 : if (ierr /= 0) then
1038 0 : call mesa_error(__FILE__,__LINE__,'init_element_atomic_weights')
1039 0 : return
1040 : end if
1041 3600 : element_atomic_weight(i) = chem_isos% W(isotope_index)
1042 : end do
1043 : end subroutine init_element_atomic_weights
1044 :
1045 :
1046 16 : subroutine init_AG89_data ! fraction by mass of total Z
1047 : ! anders & grevesse 1989, paper not available on ADS
1048 : integer :: i
1049 16 : real(dp) :: z_sum
1050 :
1051 16 : AG89_element_zfrac(:) = 0d0
1052 :
1053 16 : AG89_element_zfrac(e_li) = 5.2663415161043120D-07
1054 16 : AG89_element_zfrac(e_be) = 8.7533846996258026D-09
1055 16 : AG89_element_zfrac(e_b) = 3.0535981584981396D-07
1056 16 : AG89_element_zfrac(e_c) = 1.6164192224602159D-01
1057 16 : AG89_element_zfrac(e_n) = 5.8425953868553440D-02
1058 16 : AG89_element_zfrac(e_o) = 5.0655547566525427D-01
1059 16 : AG89_element_zfrac(e_f) = 2.1339634220190109D-05
1060 16 : AG89_element_zfrac(e_ne) = 9.2345201069972446D-02
1061 16 : AG89_element_zfrac(e_na) = 1.7588936076737712D-03
1062 16 : AG89_element_zfrac(e_mg) = 3.4766459385626718D-02
1063 16 : AG89_element_zfrac(e_al) = 3.0576538214253392D-03
1064 16 : AG89_element_zfrac(e_si) = 3.7438035164761560D-02
1065 16 : AG89_element_zfrac(e_p) = 4.2953684074804963D-04
1066 16 : AG89_element_zfrac(e_s) = 2.2000396518735602D-02
1067 16 : AG89_element_zfrac(e_cl) = 1.7836963566662123D-04
1068 16 : AG89_element_zfrac(e_ar) = 4.8868631434862610D-03
1069 16 : AG89_element_zfrac(e_k) = 1.8289986382724957D-04
1070 16 : AG89_element_zfrac(e_ca) = 3.2987013574392004D-03
1071 16 : AG89_element_zfrac(e_sc) = 2.0504272999081344D-06
1072 16 : AG89_element_zfrac(e_ti) = 1.5319766333496369D-04
1073 16 : AG89_element_zfrac(e_v) = 1.9889381301661414D-05
1074 16 : AG89_element_zfrac(e_cr) = 9.3528485499287748D-04
1075 16 : AG89_element_zfrac(e_mn) = 6.9978620325668458D-04
1076 16 : AG89_element_zfrac(e_fe) = 6.7005139944813982D-02
1077 16 : AG89_element_zfrac(e_co) = 1.7686377328884702D-04
1078 16 : AG89_element_zfrac(e_ni) = 3.8650578305194534D-03
1079 16 : AG89_element_zfrac(e_cu) = 4.4243068859971588D-05
1080 16 : AG89_element_zfrac(e_zn) = 1.0994428157112287D-04
1081 :
1082 1808 : z_sum = sum(AG89_element_zfrac(:))
1083 464 : do i = e_li, e_zn
1084 464 : AG89_element_zfrac(i) = AG89_element_zfrac(i) / z_sum
1085 : end do
1086 :
1087 16 : end subroutine init_AG89_data
1088 :
1089 :
1090 16 : subroutine init_GN93_data ! fraction by mass of total Z
1091 : ! Grevesse and Noels 1993
1092 : integer :: i
1093 16 : real(dp) :: z_sum
1094 : include 'formats'
1095 :
1096 1808 : GN93_element_zfrac(:) = -20.0d0
1097 :
1098 : !GN93_element_zfrac(e_H)=12.00d0
1099 : !GN93_element_zfrac(e_He)=10.99d0
1100 16 : GN93_element_zfrac(e_Li)=3.31d0 !meteor
1101 16 : GN93_element_zfrac(e_Be)=1.42d0 !meteor
1102 16 : GN93_element_zfrac(e_B )=2.79d0 !meteor
1103 16 : GN93_element_zfrac(e_C )=8.55d0
1104 16 : GN93_element_zfrac(e_N )=7.97d0
1105 16 : GN93_element_zfrac(e_O )=8.87d0
1106 16 : GN93_element_zfrac(e_F )=4.48d0
1107 16 : GN93_element_zfrac(e_Ne)=8.08d0
1108 16 : GN93_element_zfrac(e_Na)=6.33d0
1109 16 : GN93_element_zfrac(e_Mg)=7.58d0
1110 16 : GN93_element_zfrac(e_Al)=6.47d0
1111 16 : GN93_element_zfrac(e_Si)=7.55d0
1112 16 : GN93_element_zfrac(e_P )=5.45d0
1113 16 : GN93_element_zfrac(e_S )=7.20d0
1114 16 : GN93_element_zfrac(e_Cl)=5.28d0
1115 16 : GN93_element_zfrac(e_Ar)=6.52d0
1116 16 : GN93_element_zfrac(e_K )=5.12d0
1117 16 : GN93_element_zfrac(e_Ca)=6.36d0
1118 16 : GN93_element_zfrac(e_Sc)=3.17d0
1119 16 : GN93_element_zfrac(e_Ti)=5.02d0
1120 16 : GN93_element_zfrac(e_V )=4.00d0
1121 16 : GN93_element_zfrac(e_Cr)=5.67d0
1122 16 : GN93_element_zfrac(e_Mn)=5.39d0
1123 16 : GN93_element_zfrac(e_Fe)=7.50d0
1124 16 : GN93_element_zfrac(e_Co)=4.92d0
1125 16 : GN93_element_zfrac(e_Ni)=6.25d0
1126 16 : GN93_element_zfrac(e_Cu)=4.21d0
1127 16 : GN93_element_zfrac(e_Zn)=4.60d0
1128 16 : GN93_element_zfrac(e_Ga)=3.13d0
1129 16 : GN93_element_zfrac(e_Ge)=3.41d0
1130 16 : GN93_element_zfrac(e_As)=2.37d0
1131 16 : GN93_element_zfrac(e_Se)=3.38d0
1132 16 : GN93_element_zfrac(e_Br)=2.63d0
1133 16 : GN93_element_zfrac(e_Kr)=3.23d0
1134 16 : GN93_element_zfrac(e_Rb)=2.41d0
1135 16 : GN93_element_zfrac(e_Sr)=2.97d0
1136 16 : GN93_element_zfrac(e_Y )=2.24d0
1137 16 : GN93_element_zfrac(e_Zr)=2.60d0
1138 16 : GN93_element_zfrac(e_Nb)=1.42d0
1139 16 : GN93_element_zfrac(e_Mo)=1.92d0
1140 16 : GN93_element_zfrac(e_Ru)=1.84d0
1141 16 : GN93_element_zfrac(e_Rh)=1.12d0
1142 16 : GN93_element_zfrac(e_Pd)=1.69d0
1143 16 : GN93_element_zfrac(e_Ag)=1.24d0
1144 16 : GN93_element_zfrac(e_Cd)=1.77d0
1145 16 : GN93_element_zfrac(e_In)=0.82d0
1146 16 : GN93_element_zfrac(e_Sn)=2.14d0
1147 16 : GN93_element_zfrac(e_Sb)=1.03d0
1148 16 : GN93_element_zfrac(e_Te)=2.24d0
1149 16 : GN93_element_zfrac(e_I )=1.51d0
1150 16 : GN93_element_zfrac(e_Xe)=2.23d0
1151 16 : GN93_element_zfrac(e_Cs)=1.13d0
1152 16 : GN93_element_zfrac(e_Ba)=2.13d0
1153 16 : GN93_element_zfrac(e_La)=1.17d0
1154 16 : GN93_element_zfrac(e_Ce)=1.58d0
1155 16 : GN93_element_zfrac(e_Pr)=0.71d0
1156 16 : GN93_element_zfrac(e_Nd)=1.50d0
1157 16 : GN93_element_zfrac(e_Sm)=1.01d0
1158 16 : GN93_element_zfrac(e_Eu)=0.51d0
1159 16 : GN93_element_zfrac(e_Gd)=1.12d0
1160 16 : GN93_element_zfrac(e_Tb)=0.35d0
1161 16 : GN93_element_zfrac(e_Dy)=1.14d0
1162 16 : GN93_element_zfrac(e_Ho)=0.51d0
1163 16 : GN93_element_zfrac(e_Er)=0.93d0
1164 16 : GN93_element_zfrac(e_Tm)=0.15d0
1165 16 : GN93_element_zfrac(e_Yb)=1.08d0
1166 16 : GN93_element_zfrac(e_Lu)=0.13d0
1167 16 : GN93_element_zfrac(e_Hf)=0.88d0
1168 16 : GN93_element_zfrac(e_Ta)=-0.13d0
1169 16 : GN93_element_zfrac(e_W )=0.69d0
1170 16 : GN93_element_zfrac(e_Re)=0.28d0
1171 16 : GN93_element_zfrac(e_Os)=1.45d0
1172 16 : GN93_element_zfrac(e_Ir)=1.37d0
1173 16 : GN93_element_zfrac(e_Pt)=1.69d0
1174 16 : GN93_element_zfrac(e_Au)=0.87d0
1175 16 : GN93_element_zfrac(e_Hg)=1.17d0
1176 16 : GN93_element_zfrac(e_Tl)=0.83d0
1177 16 : GN93_element_zfrac(e_Pb)=2.06d0
1178 16 : GN93_element_zfrac(e_Bi)=0.71d0
1179 16 : GN93_element_zfrac(e_Th)=0.09d0
1180 16 : GN93_element_zfrac(e_U)=-0.50d0
1181 :
1182 : ! convert to fraction of Z by mass
1183 16 : z_sum = 0
1184 1456 : do i = e_li, e_u
1185 : GN93_element_zfrac(i) = &
1186 1440 : exp10(GN93_element_zfrac(i))*element_atomic_weight(i)
1187 1456 : z_sum = z_sum + GN93_element_zfrac(i)
1188 : end do
1189 1456 : do i = e_li, e_u
1190 1456 : GN93_element_zfrac(i) = GN93_element_zfrac(i) / z_sum
1191 : end do
1192 :
1193 16 : end subroutine init_GN93_data
1194 :
1195 :
1196 16 : subroutine init_GS98_data ! fraction by mass of total Z
1197 : ! Grevesse and Sauval 1998, Table 1
1198 : integer :: i
1199 16 : real(dp) :: z_sum
1200 : include 'formats'
1201 :
1202 1808 : GS98_element_zfrac(:) = -20.0d0
1203 :
1204 : !GS98_element_zfrac(e_H)=12.00d0
1205 : !GS98_element_zfrac(e_He)=10.93d0
1206 16 : GS98_element_zfrac(e_Li)=3.31d0 !meteor
1207 16 : GS98_element_zfrac(e_Be)=1.42d0 !meteor
1208 16 : GS98_element_zfrac(e_B)=2.79d0 !meteor
1209 16 : GS98_element_zfrac(e_C)=8.52d0
1210 16 : GS98_element_zfrac(e_N)=7.92d0
1211 16 : GS98_element_zfrac(e_O)=8.83d0
1212 16 : GS98_element_zfrac(e_F)=4.48d0
1213 16 : GS98_element_zfrac(e_Ne)=8.08d0
1214 16 : GS98_element_zfrac(e_Na)=6.32d0
1215 16 : GS98_element_zfrac(e_Mg)=7.58d0
1216 16 : GS98_element_zfrac(e_Al)=6.49d0
1217 16 : GS98_element_zfrac(e_Si)=7.56d0
1218 16 : GS98_element_zfrac(e_P)=5.56d0
1219 16 : GS98_element_zfrac(e_S)=7.20d0
1220 16 : GS98_element_zfrac(e_Cl)=5.28d0
1221 16 : GS98_element_zfrac(e_Ar)=6.40d0
1222 16 : GS98_element_zfrac(e_K)=5.13d0
1223 16 : GS98_element_zfrac(e_Ca)=6.35d0
1224 16 : GS98_element_zfrac(e_Sc)=3.10d0
1225 16 : GS98_element_zfrac(e_Ti)=4.94d0
1226 16 : GS98_element_zfrac(e_V)=4.02d0
1227 16 : GS98_element_zfrac(e_Cr)=5.69d0
1228 16 : GS98_element_zfrac(e_Mn)=5.53d0
1229 16 : GS98_element_zfrac(e_Fe)=7.50d0
1230 16 : GS98_element_zfrac(e_Co)=4.91d0
1231 16 : GS98_element_zfrac(e_Ni)=6.25d0
1232 16 : GS98_element_zfrac(e_Cu)=4.29d0
1233 16 : GS98_element_zfrac(e_Zn)=4.67d0
1234 16 : GS98_element_zfrac(e_Ga)=3.13d0
1235 16 : GS98_element_zfrac(e_Ge)=3.63d0
1236 16 : GS98_element_zfrac(e_As)=2.37d0
1237 16 : GS98_element_zfrac(e_Se)=3.41d0
1238 16 : GS98_element_zfrac(e_Br)=2.63d0
1239 16 : GS98_element_zfrac(e_Kr)=3.31d0
1240 16 : GS98_element_zfrac(e_Rb)=2.41d0
1241 16 : GS98_element_zfrac(e_Sr)=2.92d0
1242 16 : GS98_element_zfrac(e_Y)=2.23d0
1243 16 : GS98_element_zfrac(e_Zr)=2.61d0
1244 16 : GS98_element_zfrac(e_Nb)=1.40d0
1245 16 : GS98_element_zfrac(e_Mo)=1.97d0
1246 16 : GS98_element_zfrac(e_Ru)=1.83d0
1247 16 : GS98_element_zfrac(e_Rh)=1.10d0
1248 16 : GS98_element_zfrac(e_Pd)=1.70d0
1249 16 : GS98_element_zfrac(e_Ag)=1.24d0
1250 16 : GS98_element_zfrac(e_Cd)=1.76d0
1251 16 : GS98_element_zfrac(e_In)=0.82d0
1252 16 : GS98_element_zfrac(e_Sn)=2.14d0
1253 16 : GS98_element_zfrac(e_Sb)=1.03d0
1254 16 : GS98_element_zfrac(e_Te)=2.24d0
1255 16 : GS98_element_zfrac(e_I)=1.51d0
1256 16 : GS98_element_zfrac(e_Xe)=2.17d0
1257 16 : GS98_element_zfrac(e_Cs)=1.13d0
1258 16 : GS98_element_zfrac(e_Ba)=2.22d0
1259 16 : GS98_element_zfrac(e_La)=1.22d0
1260 16 : GS98_element_zfrac(e_Ce)=1.63d0
1261 16 : GS98_element_zfrac(e_Pr)=0.80d0
1262 16 : GS98_element_zfrac(e_Nd)=1.49d0
1263 16 : GS98_element_zfrac(e_Sm)=0.98d0
1264 16 : GS98_element_zfrac(e_Eu)=0.55d0
1265 16 : GS98_element_zfrac(e_Gd)=1.09d0
1266 16 : GS98_element_zfrac(e_Tb)=0.35d0
1267 16 : GS98_element_zfrac(e_Dy)=1.17d0
1268 16 : GS98_element_zfrac(e_Ho)=0.51d0
1269 16 : GS98_element_zfrac(e_Er)=0.97d0
1270 16 : GS98_element_zfrac(e_Tm)=0.15d0
1271 16 : GS98_element_zfrac(e_Yb)=0.96d0
1272 16 : GS98_element_zfrac(e_Lu)=0.13d0
1273 16 : GS98_element_zfrac(e_Hf)=0.75d0
1274 16 : GS98_element_zfrac(e_Ta)=-0.13d0
1275 16 : GS98_element_zfrac(e_W)=0.69d0
1276 16 : GS98_element_zfrac(e_Re)=0.28d0
1277 16 : GS98_element_zfrac(e_Os)=1.39d0
1278 16 : GS98_element_zfrac(e_Ir)=1.37d0
1279 16 : GS98_element_zfrac(e_Pt)=1.69d0
1280 16 : GS98_element_zfrac(e_Au)=0.85d0
1281 16 : GS98_element_zfrac(e_Hg)=1.13d0
1282 16 : GS98_element_zfrac(e_Tl)=0.83d0
1283 16 : GS98_element_zfrac(e_Pb)=2.06d0
1284 16 : GS98_element_zfrac(e_Bi)=0.71d0
1285 16 : GS98_element_zfrac(e_Th)=0.09d0
1286 16 : GS98_element_zfrac(e_U)=-0.50d0
1287 :
1288 : ! convert to fraction of Z by mass
1289 16 : z_sum = 0d0
1290 1456 : do i = e_li, e_u
1291 : GS98_element_zfrac(i) = &
1292 1440 : exp10(GS98_element_zfrac(i))*element_atomic_weight(i)
1293 1456 : z_sum = z_sum + GS98_element_zfrac(i)
1294 : end do
1295 1456 : do i = e_li, e_u
1296 1456 : GS98_element_zfrac(i) = GS98_element_zfrac(i) / z_sum
1297 : end do
1298 :
1299 16 : end subroutine init_GS98_data
1300 :
1301 16 : subroutine init_L03_data ! fraction by mass of total Z
1302 : ! Lodders 2003, ApJ, Table 1 recommended abundance
1303 : integer :: i
1304 16 : real(dp) :: z_sum
1305 : include 'formats'
1306 :
1307 1808 : L03_element_zfrac(:) = -20.0d0
1308 :
1309 : !L03_element_zfrac(e_H)=12d0
1310 : !L03_element_zfrac(e_He)=10.89d0
1311 16 : L03_element_zfrac(e_Li)=3.28d0
1312 16 : L03_element_zfrac(e_Be)=1.41d0
1313 16 : L03_element_zfrac(e_B)=2.78d0
1314 16 : L03_element_zfrac(e_C)=8.39d0
1315 16 : L03_element_zfrac(e_N)=7.83d0
1316 16 : L03_element_zfrac(e_O)=8.69d0
1317 16 : L03_element_zfrac(e_F)=4.46d0
1318 16 : L03_element_zfrac(e_Ne)=7.87d0
1319 16 : L03_element_zfrac(e_Na)=6.30d0
1320 16 : L03_element_zfrac(e_Mg)=7.55d0
1321 16 : L03_element_zfrac(e_Al)=6.46d0
1322 16 : L03_element_zfrac(e_Si)=7.54d0
1323 16 : L03_element_zfrac(e_P)=5.46d0
1324 16 : L03_element_zfrac(e_S)=7.19d0
1325 16 : L03_element_zfrac(e_Cl)=5.26d0
1326 16 : L03_element_zfrac(e_Ar)=6.55d0
1327 16 : L03_element_zfrac(e_K)=5.11d0
1328 16 : L03_element_zfrac(e_Ca)=6.34d0
1329 16 : L03_element_zfrac(e_Sc)=3.07d0
1330 16 : L03_element_zfrac(e_Ti)=4.92d0
1331 16 : L03_element_zfrac(e_V)=4.00d0
1332 16 : L03_element_zfrac(e_Cr)=5.65d0
1333 16 : L03_element_zfrac(e_Mn)=5.50d0
1334 16 : L03_element_zfrac(e_Fe)=7.47d0
1335 16 : L03_element_zfrac(e_Co)=4.91d0
1336 16 : L03_element_zfrac(e_Ni)=6.22d0
1337 16 : L03_element_zfrac(e_Cu)=4.26d0
1338 16 : L03_element_zfrac(e_Zn)=4.63d0
1339 16 : L03_element_zfrac(e_Ga)=3.10d0
1340 16 : L03_element_zfrac(e_Ge)=3.62d0
1341 16 : L03_element_zfrac(e_As)=2.32d0
1342 16 : L03_element_zfrac(e_Se)=3.36d0
1343 16 : L03_element_zfrac(e_Br)=2.59d0
1344 16 : L03_element_zfrac(e_Kr)=3.28d0
1345 16 : L03_element_zfrac(e_Rb)=2.36d0
1346 16 : L03_element_zfrac(e_Sr)=2.91d0
1347 16 : L03_element_zfrac(e_Y)=2.20d0
1348 16 : L03_element_zfrac(e_Zr)=2.60d0
1349 16 : L03_element_zfrac(e_Nb)=1.42d0
1350 16 : L03_element_zfrac(e_Mo)=1.96d0
1351 16 : L03_element_zfrac(e_Ru)=1.82d0
1352 16 : L03_element_zfrac(e_Rh)=1.11d0
1353 16 : L03_element_zfrac(e_Pd)=1.70d0
1354 16 : L03_element_zfrac(e_Ag)=1.23d0
1355 16 : L03_element_zfrac(e_Cd)=1.74d0
1356 16 : L03_element_zfrac(e_In)=0.80d0
1357 16 : L03_element_zfrac(e_Sn)=2.11d0
1358 16 : L03_element_zfrac(e_Sb)=1.06d0
1359 16 : L03_element_zfrac(e_Te)=2.22d0
1360 16 : L03_element_zfrac(e_I)=1.54d0
1361 16 : L03_element_zfrac(e_Xe)=2.27d0
1362 16 : L03_element_zfrac(e_Cs)=1.10d0
1363 16 : L03_element_zfrac(e_Ba)=2.18d0
1364 16 : L03_element_zfrac(e_La)=1.18d0
1365 16 : L03_element_zfrac(e_Ce)=1.61d0
1366 16 : L03_element_zfrac(e_Pr)=0.78d0
1367 16 : L03_element_zfrac(e_Nd)=1.46d0
1368 16 : L03_element_zfrac(e_Sm)=0.95d0
1369 16 : L03_element_zfrac(e_Eu)=0.52d0
1370 16 : L03_element_zfrac(e_Gd)=1.06d0
1371 16 : L03_element_zfrac(e_Tb)=0.31d0
1372 16 : L03_element_zfrac(e_Dy)=1.13d0
1373 16 : L03_element_zfrac(e_Ho)=0.49d0
1374 16 : L03_element_zfrac(e_Er)=0.95d0
1375 16 : L03_element_zfrac(e_Tm)=0.11d0
1376 16 : L03_element_zfrac(e_Yb)=0.94d0
1377 16 : L03_element_zfrac(e_Lu)=0.09d0
1378 16 : L03_element_zfrac(e_Hf)=0.77d0
1379 16 : L03_element_zfrac(e_Ta)=-0.14d0
1380 16 : L03_element_zfrac(e_W)=0.65d0
1381 16 : L03_element_zfrac(e_Re)=0.26d0
1382 16 : L03_element_zfrac(e_Os)=1.37d0
1383 16 : L03_element_zfrac(e_Ir)=1.35d0
1384 16 : L03_element_zfrac(e_Pt)=1.67d0
1385 16 : L03_element_zfrac(e_Au)=0.83d0
1386 16 : L03_element_zfrac(e_Hg)=1.16d0
1387 16 : L03_element_zfrac(e_Tl)=0.81d0
1388 16 : L03_element_zfrac(e_Pb)=2.05d0
1389 16 : L03_element_zfrac(e_Bi)=0.68d0
1390 16 : L03_element_zfrac(e_Th)=0.09d0
1391 16 : L03_element_zfrac(e_U)=-0.49d0
1392 :
1393 : ! convert to fraction of Z by mass
1394 16 : z_sum = 0d0
1395 1456 : do i = e_li, e_u
1396 : L03_element_zfrac(i) = &
1397 1440 : exp10(L03_element_zfrac(i))*element_atomic_weight(i)
1398 1456 : z_sum = z_sum + L03_element_zfrac(i)
1399 : end do
1400 1456 : do i = e_li, e_u
1401 1456 : L03_element_zfrac(i) = L03_element_zfrac(i) / z_sum
1402 : end do
1403 :
1404 16 : end subroutine init_L03_data
1405 :
1406 :
1407 16 : subroutine init_AGS05_data ! fraction by mass of total Z
1408 : ! Asplund, Grevesse and Sauval 2005
1409 : integer :: i
1410 16 : real(dp) :: z_sum
1411 : include 'formats'
1412 :
1413 1808 : AGS05_element_zfrac(:) = -20.0d0
1414 :
1415 : ! first store log abundances from the paper (photosphere unless otherwise noted)
1416 : ! relative to log abundance of H = 12.00d0
1417 16 : AGS05_element_zfrac(e_li) = 3.25d0 !meteor
1418 16 : AGS05_element_zfrac(e_be) = 1.38d0
1419 16 : AGS05_element_zfrac(e_b ) = 2.70d0
1420 16 : AGS05_element_zfrac(e_c ) = 8.39d0
1421 16 : AGS05_element_zfrac(e_n ) = 7.78d0
1422 16 : AGS05_element_zfrac(e_o ) = 8.66d0
1423 16 : AGS05_element_zfrac(e_f ) = 4.56d0
1424 16 : AGS05_element_zfrac(e_ne) = 7.84d0 !indirect
1425 16 : AGS05_element_zfrac(e_na) = 6.17d0
1426 16 : AGS05_element_zfrac(e_mg) = 7.53d0
1427 16 : AGS05_element_zfrac(e_al) = 6.37d0
1428 16 : AGS05_element_zfrac(e_si) = 7.51d0
1429 16 : AGS05_element_zfrac(e_p ) = 5.36d0
1430 16 : AGS05_element_zfrac(e_s ) = 7.14d0
1431 16 : AGS05_element_zfrac(e_cl) = 5.50d0
1432 16 : AGS05_element_zfrac(e_ar) = 6.18d0 !indirect
1433 16 : AGS05_element_zfrac(e_k ) = 5.08d0
1434 16 : AGS05_element_zfrac(e_ca) = 6.31d0
1435 16 : AGS05_element_zfrac(e_sc) = 3.05d0
1436 16 : AGS05_element_zfrac(e_ti) = 4.90d0
1437 16 : AGS05_element_zfrac(e_v ) = 4.00d0
1438 16 : AGS05_element_zfrac(e_cr) = 5.64d0
1439 16 : AGS05_element_zfrac(e_mn) = 5.39d0
1440 16 : AGS05_element_zfrac(e_fe) = 7.45d0
1441 16 : AGS05_element_zfrac(e_co) = 4.92d0
1442 16 : AGS05_element_zfrac(e_ni) = 6.23d0
1443 16 : AGS05_element_zfrac(e_cu) = 4.21d0
1444 16 : AGS05_element_zfrac(e_zn) = 4.60d0
1445 16 : AGS05_element_zfrac(e_ga) = 2.88d0
1446 16 : AGS05_element_zfrac(e_ge) = 3.58d0
1447 16 : AGS05_element_zfrac(e_as) = 2.29d0 !meteor
1448 16 : AGS05_element_zfrac(e_se) = 3.33d0 !meteor
1449 16 : AGS05_element_zfrac(e_br) = 2.56d0 !meteor
1450 16 : AGS05_element_zfrac(e_kr) = 3.28d0 !indirect
1451 16 : AGS05_element_zfrac(e_rb) = 2.60d0
1452 16 : AGS05_element_zfrac(e_sr) = 2.92d0
1453 16 : AGS05_element_zfrac(e_y ) = 2.21d0
1454 16 : AGS05_element_zfrac(e_zr) = 2.59d0
1455 16 : AGS05_element_zfrac(e_nb) = 1.42d0
1456 16 : AGS05_element_zfrac(e_mo) = 1.92d0
1457 16 : AGS05_element_zfrac(e_Ru) = 1.84d0
1458 16 : AGS05_element_zfrac(e_Rh) = 1.12d0
1459 16 : AGS05_element_zfrac(e_Pd) = 1.69d0
1460 16 : AGS05_element_zfrac(e_Ag) = 0.94d0
1461 16 : AGS05_element_zfrac(e_Cd) = 1.77d0
1462 16 : AGS05_element_zfrac(e_In) = 1.60d0
1463 16 : AGS05_element_zfrac(e_Sn) = 2.00d0
1464 16 : AGS05_element_zfrac(e_Sb) = 1.00d0
1465 16 : AGS05_element_zfrac(e_Te) = 2.19d0 !meteor
1466 16 : AGS05_element_zfrac(e_I ) = 1.51d0 !meteor
1467 16 : AGS05_element_zfrac(e_Xe) = 2.27d0 !indirect
1468 16 : AGS05_element_zfrac(e_Cs) = 1.07d0 !meteor
1469 16 : AGS05_element_zfrac(e_Ba) = 2.17d0
1470 16 : AGS05_element_zfrac(e_La) = 1.13d0
1471 16 : AGS05_element_zfrac(e_Ce) = 1.58d0
1472 16 : AGS05_element_zfrac(e_Pr) = 0.71d0
1473 16 : AGS05_element_zfrac(e_Nd) = 1.45d0
1474 16 : AGS05_element_zfrac(e_Sm) = 1.01d0
1475 16 : AGS05_element_zfrac(e_Eu) = 0.52d0
1476 16 : AGS05_element_zfrac(e_Gd) = 1.12d0
1477 16 : AGS05_element_zfrac(e_Tb) = 0.28d0
1478 16 : AGS05_element_zfrac(e_Dy) = 1.14d0
1479 16 : AGS05_element_zfrac(e_Ho) = 0.51d0
1480 16 : AGS05_element_zfrac(e_Er) = 0.93d0
1481 16 : AGS05_element_zfrac(e_Tm) = 0.00d0
1482 16 : AGS05_element_zfrac(e_Yb) = 1.08d0
1483 16 : AGS05_element_zfrac(e_Lu) = 0.06d0
1484 16 : AGS05_element_zfrac(e_Hf) = 0.88d0
1485 16 : AGS05_element_zfrac(e_Ta) = -0.17d0 !meteor
1486 16 : AGS05_element_zfrac(e_W ) = 1.11d0
1487 16 : AGS05_element_zfrac(e_Re) = 0.23d0 !meteor
1488 16 : AGS05_element_zfrac(e_Os) = 1.45d0
1489 16 : AGS05_element_zfrac(e_Ir) = 1.38d0
1490 16 : AGS05_element_zfrac(e_Pt) = 1.64d0 !meteor
1491 16 : AGS05_element_zfrac(e_Au) = 1.01d0
1492 16 : AGS05_element_zfrac(e_Hg) = 1.13d0 !meteor
1493 16 : AGS05_element_zfrac(e_Tl) = 0.90d0
1494 16 : AGS05_element_zfrac(e_Pb) = 2.00d0
1495 16 : AGS05_element_zfrac(e_Bi) = 0.65d0 !meteor
1496 16 : AGS05_element_zfrac(e_Th) = 0.06d0 !meteor
1497 16 : AGS05_element_zfrac(e_U) = -0.52d0
1498 :
1499 : ! convert to fraction of Z by mass
1500 16 : z_sum = 0d0
1501 1456 : do i = e_li, e_u
1502 : AGS05_element_zfrac(i) = &
1503 1440 : exp10(AGS05_element_zfrac(i))*element_atomic_weight(i)
1504 1456 : z_sum = z_sum + AGS05_element_zfrac(i)
1505 : end do
1506 1456 : do i = e_li, e_u
1507 1456 : AGS05_element_zfrac(i) = AGS05_element_zfrac(i) / z_sum
1508 : end do
1509 :
1510 16 : end subroutine init_AGS05_data
1511 :
1512 :
1513 16 : subroutine init_AGSS09_data ! fraction by mass of total Z
1514 : ! Asplund, Grevesse, Sauval, and Scott 2009 abundances
1515 : ! Annu. Rev. Astron. Astrophys. 2009. 47:481–522
1516 : integer :: i
1517 16 : real(dp) :: z_sum
1518 : include 'formats'
1519 :
1520 1808 : AGSS09_element_zfrac(:) = -20.0d0
1521 :
1522 : ! first store log abundances from the paper
1523 16 : AGSS09_element_zfrac(e_li) = 3.26d0
1524 16 : AGSS09_element_zfrac(e_be) = 1.38d0
1525 16 : AGSS09_element_zfrac(e_b ) = 2.70d0
1526 16 : AGSS09_element_zfrac(e_c ) = 8.43d0
1527 16 : AGSS09_element_zfrac(e_n ) = 7.83d0
1528 16 : AGSS09_element_zfrac(e_o ) = 8.69d0
1529 16 : AGSS09_element_zfrac(e_f ) = 4.56d0
1530 16 : AGSS09_element_zfrac(e_ne) = 7.93d0
1531 16 : AGSS09_element_zfrac(e_na) = 6.24d0
1532 16 : AGSS09_element_zfrac(e_mg) = 7.60d0
1533 16 : AGSS09_element_zfrac(e_al) = 6.45d0
1534 16 : AGSS09_element_zfrac(e_si) = 7.51d0
1535 16 : AGSS09_element_zfrac(e_p ) = 5.41d0
1536 16 : AGSS09_element_zfrac(e_s ) = 7.12d0
1537 16 : AGSS09_element_zfrac(e_cl) = 5.50d0
1538 16 : AGSS09_element_zfrac(e_ar) = 6.40d0
1539 16 : AGSS09_element_zfrac(e_k ) = 5.03d0
1540 16 : AGSS09_element_zfrac(e_ca) = 6.34d0
1541 16 : AGSS09_element_zfrac(e_sc) = 3.15d0
1542 16 : AGSS09_element_zfrac(e_ti) = 4.95d0
1543 16 : AGSS09_element_zfrac(e_v ) = 3.93d0
1544 16 : AGSS09_element_zfrac(e_cr) = 5.64d0
1545 16 : AGSS09_element_zfrac(e_mn) = 5.43d0
1546 16 : AGSS09_element_zfrac(e_fe) = 7.50d0
1547 16 : AGSS09_element_zfrac(e_co) = 4.99d0
1548 16 : AGSS09_element_zfrac(e_ni) = 6.22d0
1549 16 : AGSS09_element_zfrac(e_cu) = 4.19d0
1550 16 : AGSS09_element_zfrac(e_zn) = 4.56d0
1551 16 : AGSS09_element_zfrac(e_ga) = 3.04d0
1552 16 : AGSS09_element_zfrac(e_ge) = 3.65d0
1553 16 : AGSS09_element_zfrac(e_as) = 2.30d0 !meteor
1554 16 : AGSS09_element_zfrac(e_se) = 3.34d0 !meteor
1555 16 : AGSS09_element_zfrac(e_br) = 2.54d0 !meteor
1556 16 : AGSS09_element_zfrac(e_kr) = 3.25d0 !indirect
1557 16 : AGSS09_element_zfrac(e_rb) = 2.52d0
1558 16 : AGSS09_element_zfrac(e_sr) = 2.87d0
1559 16 : AGSS09_element_zfrac(e_y ) = 2.21d0
1560 16 : AGSS09_element_zfrac(e_zr) = 2.58d0
1561 16 : AGSS09_element_zfrac(e_nb) = 1.46d0
1562 16 : AGSS09_element_zfrac(e_mo) = 1.88d0
1563 16 : AGSS09_element_zfrac(e_Ru) = 1.75d0
1564 16 : AGSS09_element_zfrac(e_Rh) = 0.91d0
1565 16 : AGSS09_element_zfrac(e_Pd) = 1.57d0
1566 16 : AGSS09_element_zfrac(e_Ag) = 0.94d0
1567 16 : AGSS09_element_zfrac(e_Cd) = 1.71d0
1568 16 : AGSS09_element_zfrac(e_In) = 0.80d0
1569 16 : AGSS09_element_zfrac(e_Sn) = 2.04d0
1570 16 : AGSS09_element_zfrac(e_Sb) = 1.01d0
1571 16 : AGSS09_element_zfrac(e_Te) = 2.18d0
1572 16 : AGSS09_element_zfrac(e_I ) = 1.55d0
1573 16 : AGSS09_element_zfrac(e_Xe) = 2.24d0
1574 16 : AGSS09_element_zfrac(e_Cs) = 1.08d0
1575 16 : AGSS09_element_zfrac(e_Ba) = 2.18d0
1576 16 : AGSS09_element_zfrac(e_La) = 1.10d0
1577 16 : AGSS09_element_zfrac(e_Ce) = 1.58d0
1578 16 : AGSS09_element_zfrac(e_Pr) = 0.72d0
1579 16 : AGSS09_element_zfrac(e_Nd) = 1.42d0
1580 16 : AGSS09_element_zfrac(e_Sm) = 0.96d0
1581 16 : AGSS09_element_zfrac(e_Eu) = 0.52d0
1582 16 : AGSS09_element_zfrac(e_Gd) = 1.07d0
1583 16 : AGSS09_element_zfrac(e_Tb) = 0.30d0
1584 16 : AGSS09_element_zfrac(e_Dy) = 1.10d0
1585 16 : AGSS09_element_zfrac(e_Ho) = 0.48d0
1586 16 : AGSS09_element_zfrac(e_Er) = 0.92d0
1587 16 : AGSS09_element_zfrac(e_Tm) = 0.10d0
1588 16 : AGSS09_element_zfrac(e_Yb) = 0.84d0
1589 16 : AGSS09_element_zfrac(e_Lu) = 0.10d0
1590 16 : AGSS09_element_zfrac(e_Hf) = 0.85d0
1591 16 : AGSS09_element_zfrac(e_Ta) = -0.12d0
1592 16 : AGSS09_element_zfrac(e_W ) = 0.85d0
1593 16 : AGSS09_element_zfrac(e_Re) = 0.26d0
1594 16 : AGSS09_element_zfrac(e_Os) = 1.40d0
1595 16 : AGSS09_element_zfrac(e_Ir) = 1.38d0
1596 16 : AGSS09_element_zfrac(e_Pt) = 1.62d0
1597 16 : AGSS09_element_zfrac(e_Au) = 0.92d0
1598 16 : AGSS09_element_zfrac(e_Hg) = 1.17d0
1599 16 : AGSS09_element_zfrac(e_Tl) = 0.90d0
1600 16 : AGSS09_element_zfrac(e_Pb) = 1.75d0
1601 16 : AGSS09_element_zfrac(e_Bi) = 0.65d0
1602 16 : AGSS09_element_zfrac(e_Th) = 0.02d0
1603 16 : AGSS09_element_zfrac(e_U) = -0.54d0
1604 :
1605 : ! convert to fraction of Z by mass
1606 16 : z_sum = 0
1607 1456 : do i = e_li, e_u
1608 : AGSS09_element_zfrac(i) = &
1609 1440 : exp10(AGSS09_element_zfrac(i))*element_atomic_weight(i)
1610 1456 : z_sum = z_sum + AGSS09_element_zfrac(i)
1611 : end do
1612 1456 : do i = e_li, e_u
1613 1456 : AGSS09_element_zfrac(i) = AGSS09_element_zfrac(i) / z_sum
1614 : end do
1615 :
1616 16 : end subroutine init_AGSS09_data
1617 :
1618 :
1619 16 : subroutine init_A09_Przybilla_data ! fraction by mass of total Z
1620 : ! provided by Ehsan Moravveji, Oct 12, 2013.
1621 : ! The mass fraction is taken from Asplund et al. (2009), and modified by the
1622 : ! B-star measurement of Nieva & Przybilla 2012, A&A, 539, 143 and
1623 : ! Przybilla et al. (2013), EAS proceeding to be published
1624 : ! The modified elements are: he, c, n, o, ne, mg, al, si, s, ar, fe
1625 : integer :: i
1626 16 : real(dp) :: z_sum
1627 : include 'formats'
1628 :
1629 1808 : A09_Prz_zfrac(:) = -20.0d0
1630 :
1631 : ! A09_Prz_zfrac(e_h ) = 12.00d0
1632 : ! A09_Prz_zfrac(e_he) = 10.99d0
1633 16 : A09_Prz_zfrac(e_li) = 3.26d0
1634 16 : A09_Prz_zfrac(e_be) = 1.38d0
1635 16 : A09_Prz_zfrac(e_b ) = 2.70d0
1636 16 : A09_Prz_zfrac(e_c ) = 8.33d0
1637 16 : A09_Prz_zfrac(e_n ) = 7.79d0
1638 16 : A09_Prz_zfrac(e_o ) = 8.76d0
1639 16 : A09_Prz_zfrac(e_f ) = 4.56d0
1640 16 : A09_Prz_zfrac(e_ne) = 8.09d0
1641 16 : A09_Prz_zfrac(e_na) = 6.24d0
1642 16 : A09_Prz_zfrac(e_mg) = 7.56d0
1643 16 : A09_Prz_zfrac(e_al) = 6.30d0
1644 16 : A09_Prz_zfrac(e_si) = 7.50d0
1645 16 : A09_Prz_zfrac(e_p ) = 5.41d0
1646 16 : A09_Prz_zfrac(e_s ) = 7.14d0
1647 16 : A09_Prz_zfrac(e_cl) = 5.50d0
1648 16 : A09_Prz_zfrac(e_ar) = 6.50d0
1649 16 : A09_Prz_zfrac(e_k ) = 5.03d0
1650 16 : A09_Prz_zfrac(e_ca) = 6.34d0
1651 16 : A09_Prz_zfrac(e_sc) = 3.15d0
1652 16 : A09_Prz_zfrac(e_ti) = 4.95d0
1653 16 : A09_Prz_zfrac(e_v ) = 3.93d0
1654 16 : A09_Prz_zfrac(e_cr) = 5.64d0
1655 16 : A09_Prz_zfrac(e_mn) = 5.43d0
1656 16 : A09_Prz_zfrac(e_fe) = 7.52d0
1657 16 : A09_Prz_zfrac(e_co) = 4.99d0
1658 16 : A09_Prz_zfrac(e_ni) = 6.22d0
1659 16 : A09_Prz_zfrac(e_cu) = 4.19d0
1660 16 : A09_Prz_zfrac(e_zn) = 4.56d0
1661 16 : A09_Prz_zfrac(e_ga) = 3.04d0
1662 16 : A09_Prz_zfrac(e_ge) = 3.65d0
1663 16 : A09_Prz_zfrac(e_as) = 2.30d0
1664 16 : A09_Prz_zfrac(e_se) = 3.34d0
1665 16 : A09_Prz_zfrac(e_br) = 2.54d0
1666 16 : A09_Prz_zfrac(e_kr) = 3.25d0
1667 16 : A09_Prz_zfrac(e_rb) = 2.52d0
1668 16 : A09_Prz_zfrac(e_sr) = 2.87d0
1669 16 : A09_Prz_zfrac(e_y ) = 2.21d0
1670 16 : A09_Prz_zfrac(e_zr) = 2.58d0
1671 16 : A09_Prz_zfrac(e_nb) = 1.46d0
1672 16 : A09_Prz_zfrac(e_mo) = 1.88d0
1673 16 : A09_Prz_zfrac(e_ru) = 1.75d0
1674 16 : A09_Prz_zfrac(e_rh) = 0.91d0
1675 16 : A09_Prz_zfrac(e_pd) = 1.57d0
1676 16 : A09_Prz_zfrac(e_ag) = 0.94d0
1677 16 : A09_Prz_zfrac(e_cd) = 1.71d0
1678 16 : A09_Prz_zfrac(e_in) = 0.80d0
1679 16 : A09_Prz_zfrac(e_sn) = 2.04d0
1680 16 : A09_Prz_zfrac(e_sb) = 1.01d0
1681 16 : A09_Prz_zfrac(e_te) = 2.18d0
1682 16 : A09_Prz_zfrac(e_i ) = 1.55d0
1683 16 : A09_Prz_zfrac(e_xe) = 2.24d0
1684 16 : A09_Prz_zfrac(e_cs) = 1.08d0
1685 16 : A09_Prz_zfrac(e_ba) = 2.18d0
1686 16 : A09_Prz_zfrac(e_la) = 1.10d0
1687 16 : A09_Prz_zfrac(e_ce) = 1.58d0
1688 16 : A09_Prz_zfrac(e_pr) = 0.72d0
1689 16 : A09_Prz_zfrac(e_nd) = 1.42d0
1690 16 : A09_Prz_zfrac(e_sm) = 0.96d0
1691 16 : A09_Prz_zfrac(e_eu) = 0.52d0
1692 16 : A09_Prz_zfrac(e_gd) = 1.07d0
1693 16 : A09_Prz_zfrac(e_tb) = 0.30d0
1694 16 : A09_Prz_zfrac(e_dy) = 1.10d0
1695 16 : A09_Prz_zfrac(e_ho) = 0.48d0
1696 16 : A09_Prz_zfrac(e_er) = 0.92d0
1697 16 : A09_Prz_zfrac(e_tm) = 0.10d0
1698 16 : A09_Prz_zfrac(e_yb) = 0.84d0
1699 16 : A09_Prz_zfrac(e_lu) = 0.10d0
1700 16 : A09_Prz_zfrac(e_hf) = 0.85d0
1701 16 : A09_Prz_zfrac(e_ta) = -0.12d0
1702 16 : A09_Prz_zfrac(e_w ) = 0.85d0
1703 16 : A09_Prz_zfrac(e_re) = 0.26d0
1704 16 : A09_Prz_zfrac(e_os) = 1.40d0
1705 16 : A09_Prz_zfrac(e_ir) = 1.38d0
1706 16 : A09_Prz_zfrac(e_pt) = 1.62d0
1707 16 : A09_Prz_zfrac(e_au) = 0.92d0
1708 16 : A09_Prz_zfrac(e_hg) = 1.17d0
1709 16 : A09_Prz_zfrac(e_tl) = 0.90d0
1710 16 : A09_Prz_zfrac(e_pb) = 1.75d0
1711 16 : A09_Prz_zfrac(e_bi) = 0.65d0
1712 16 : A09_Prz_zfrac(e_th) = 0.02d0
1713 16 : A09_Prz_zfrac(e_u ) = -0.54d0
1714 :
1715 : ! convert to fraction of Z by mass
1716 16 : z_sum = 0d0
1717 1456 : do i = e_li, e_u
1718 : A09_Prz_zfrac(i) = &
1719 1440 : exp10(A09_Prz_zfrac(i))*element_atomic_weight(i)
1720 1456 : z_sum = z_sum + A09_Prz_zfrac(i)
1721 : end do
1722 1456 : do i = e_li, e_u
1723 1456 : A09_Prz_zfrac(i) = A09_Prz_zfrac(i) / z_sum
1724 : end do
1725 :
1726 16 : end subroutine init_A09_Przybilla_data
1727 :
1728 16 : subroutine init_MB22_photospheric_data ! fraction by mass of total Z
1729 : ! Ekaterina Magg et al. , A&A 661, A140 (2022) photospheric abundance.
1730 : ! supplemented with Asplund, Grevesse, Sauval, and Scott 2009 abundances
1731 : integer :: i
1732 16 : real(dp) :: z_sum
1733 : include 'formats'
1734 :
1735 1808 : MB22_photospheric_element_zfrac(:) = -20.0d0
1736 :
1737 : ! first store log abundances from the paper
1738 16 : MB22_photospheric_element_zfrac(e_li) = 3.26d0 ! AGSS09
1739 16 : MB22_photospheric_element_zfrac(e_be) = 1.38d0 ! AGSS09
1740 16 : MB22_photospheric_element_zfrac(e_b ) = 2.70d0 ! AGSS09
1741 16 : MB22_photospheric_element_zfrac(e_c ) = 8.56d0
1742 16 : MB22_photospheric_element_zfrac(e_n ) = 7.98d0
1743 16 : MB22_photospheric_element_zfrac(e_o ) = 8.77d0
1744 16 : MB22_photospheric_element_zfrac(e_f ) = 4.40d0
1745 16 : MB22_photospheric_element_zfrac(e_ne) = 8.15d0
1746 16 : MB22_photospheric_element_zfrac(e_na) = 6.29d0
1747 16 : MB22_photospheric_element_zfrac(e_mg) = 7.55d0
1748 16 : MB22_photospheric_element_zfrac(e_al) = 6.43d0
1749 16 : MB22_photospheric_element_zfrac(e_si) = 7.59d0
1750 16 : MB22_photospheric_element_zfrac(e_p ) = 5.41d0
1751 16 : MB22_photospheric_element_zfrac(e_s ) = 7.16d0
1752 16 : MB22_photospheric_element_zfrac(e_cl) = 5.25d0
1753 16 : MB22_photospheric_element_zfrac(e_ar) = 6.50d0
1754 16 : MB22_photospheric_element_zfrac(e_k ) = 5.14d0
1755 16 : MB22_photospheric_element_zfrac(e_ca) = 6.37d0
1756 16 : MB22_photospheric_element_zfrac(e_sc) = 3.07d0
1757 16 : MB22_photospheric_element_zfrac(e_ti) = 4.94d0
1758 16 : MB22_photospheric_element_zfrac(e_v ) = 3.89d0
1759 16 : MB22_photospheric_element_zfrac(e_cr) = 5.74d0
1760 16 : MB22_photospheric_element_zfrac(e_mn) = 5.52d0
1761 16 : MB22_photospheric_element_zfrac(e_fe) = 7.50d0
1762 16 : MB22_photospheric_element_zfrac(e_co) = 4.95d0
1763 16 : MB22_photospheric_element_zfrac(e_ni) = 6.24d0
1764 16 : MB22_photospheric_element_zfrac(e_cu) = 4.19d0 ! AGSS09
1765 16 : MB22_photospheric_element_zfrac(e_zn) = 4.56d0 ! AGSS09
1766 16 : MB22_photospheric_element_zfrac(e_ga) = 3.04d0 ! AGSS09
1767 16 : MB22_photospheric_element_zfrac(e_ge) = 3.65d0 ! AGSS09
1768 16 : MB22_photospheric_element_zfrac(e_as) = 2.30d0 !meteor AGSS09
1769 16 : MB22_photospheric_element_zfrac(e_se) = 3.34d0 !meteor AGSS09
1770 16 : MB22_photospheric_element_zfrac(e_br) = 2.54d0 !meteor AGSS09
1771 16 : MB22_photospheric_element_zfrac(e_kr) = 3.25d0 !indirect AGSS09
1772 16 : MB22_photospheric_element_zfrac(e_rb) = 2.52d0 ! AGSS09
1773 16 : MB22_photospheric_element_zfrac(e_sr) = 2.87d0 ! AGSS09
1774 16 : MB22_photospheric_element_zfrac(e_y ) = 2.21d0 ! AGSS09
1775 16 : MB22_photospheric_element_zfrac(e_zr) = 2.58d0 ! AGSS09
1776 16 : MB22_photospheric_element_zfrac(e_nb) = 1.46d0 ! AGSS09
1777 16 : MB22_photospheric_element_zfrac(e_mo) = 1.88d0 ! AGSS09
1778 16 : MB22_photospheric_element_zfrac(e_Ru) = 1.75d0 ! AGSS09
1779 16 : MB22_photospheric_element_zfrac(e_Rh) = 0.91d0 ! AGSS09
1780 16 : MB22_photospheric_element_zfrac(e_Pd) = 1.57d0 ! AGSS09
1781 16 : MB22_photospheric_element_zfrac(e_Ag) = 0.94d0 ! AGSS09
1782 16 : MB22_photospheric_element_zfrac(e_Cd) = 1.71d0 ! AGSS09
1783 16 : MB22_photospheric_element_zfrac(e_In) = 0.80d0 ! AGSS09
1784 16 : MB22_photospheric_element_zfrac(e_Sn) = 2.04d0 ! AGSS09
1785 16 : MB22_photospheric_element_zfrac(e_Sb) = 1.01d0 ! AGSS09
1786 16 : MB22_photospheric_element_zfrac(e_Te) = 2.18d0 ! AGSS09
1787 16 : MB22_photospheric_element_zfrac(e_I ) = 1.55d0 ! AGSS09
1788 16 : MB22_photospheric_element_zfrac(e_Xe) = 2.24d0 ! AGSS09
1789 16 : MB22_photospheric_element_zfrac(e_Cs) = 1.08d0 ! AGSS09
1790 16 : MB22_photospheric_element_zfrac(e_Ba) = 2.18d0 ! AGSS09
1791 16 : MB22_photospheric_element_zfrac(e_La) = 1.10d0 ! AGSS09
1792 16 : MB22_photospheric_element_zfrac(e_Ce) = 1.58d0 ! AGSS09
1793 16 : MB22_photospheric_element_zfrac(e_Pr) = 0.72d0 ! AGSS09
1794 16 : MB22_photospheric_element_zfrac(e_Nd) = 1.42d0 ! AGSS09
1795 16 : MB22_photospheric_element_zfrac(e_Sm) = 0.96d0 ! AGSS09
1796 16 : MB22_photospheric_element_zfrac(e_Eu) = 0.52d0 ! AGSS09
1797 16 : MB22_photospheric_element_zfrac(e_Gd) = 1.07d0 ! AGSS09
1798 16 : MB22_photospheric_element_zfrac(e_Tb) = 0.30d0 ! AGSS09
1799 16 : MB22_photospheric_element_zfrac(e_Dy) = 1.10d0 ! AGSS09
1800 16 : MB22_photospheric_element_zfrac(e_Ho) = 0.48d0 ! AGSS09
1801 16 : MB22_photospheric_element_zfrac(e_Er) = 0.92d0 ! AGSS09
1802 16 : MB22_photospheric_element_zfrac(e_Tm) = 0.10d0 ! AGSS09
1803 16 : MB22_photospheric_element_zfrac(e_Yb) = 0.84d0 ! AGSS09
1804 16 : MB22_photospheric_element_zfrac(e_Lu) = 0.10d0 ! AGSS09
1805 16 : MB22_photospheric_element_zfrac(e_Hf) = 0.85d0 ! AGSS09
1806 16 : MB22_photospheric_element_zfrac(e_Ta) = -0.12d0 ! AGSS09
1807 16 : MB22_photospheric_element_zfrac(e_W ) = 0.85d0 ! AGSS09
1808 16 : MB22_photospheric_element_zfrac(e_Re) = 0.26d0 ! AGSS09
1809 16 : MB22_photospheric_element_zfrac(e_Os) = 1.40d0 ! AGSS09
1810 16 : MB22_photospheric_element_zfrac(e_Ir) = 1.38d0 ! AGSS09
1811 16 : MB22_photospheric_element_zfrac(e_Pt) = 1.62d0 ! AGSS09
1812 16 : MB22_photospheric_element_zfrac(e_Au) = 0.92d0 ! AGSS09
1813 16 : MB22_photospheric_element_zfrac(e_Hg) = 1.17d0 ! AGSS09
1814 16 : MB22_photospheric_element_zfrac(e_Tl) = 0.90d0 ! AGSS09
1815 16 : MB22_photospheric_element_zfrac(e_Pb) = 1.75d0 ! AGSS09
1816 16 : MB22_photospheric_element_zfrac(e_Bi) = 0.65d0 ! AGSS09
1817 16 : MB22_photospheric_element_zfrac(e_Th) = 0.02d0 ! AGSS09
1818 16 : MB22_photospheric_element_zfrac(e_U) = -0.54d0 ! AGSS09
1819 :
1820 : ! convert to fraction of Z by mass
1821 16 : z_sum = 0
1822 1456 : do i = e_li, e_u
1823 : MB22_photospheric_element_zfrac(i) = &
1824 1440 : exp10(MB22_photospheric_element_zfrac(i))*element_atomic_weight(i)
1825 1456 : z_sum = z_sum + MB22_photospheric_element_zfrac(i)
1826 : end do
1827 1456 : do i = e_li, e_u
1828 1456 : MB22_photospheric_element_zfrac(i) = MB22_photospheric_element_zfrac(i) / z_sum
1829 : end do
1830 :
1831 16 : end subroutine init_MB22_photospheric_data
1832 :
1833 :
1834 16 : subroutine init_AAG21_photospheric_data ! fraction by mass of total Z
1835 : ! Asplund et al. A&A 653, A141 (2021) photospheric abundance.
1836 : ! Supplemented with meteoric values
1837 : integer :: i
1838 16 : real(dp) :: z_sum
1839 : include 'formats'
1840 :
1841 1808 : AAG21_photospheric_element_zfrac(:) = -20.0d0
1842 :
1843 : ! first store log abundances from the paper
1844 16 : AAG21_photospheric_element_zfrac(e_li) = 0.96d0 !± 0.013
1845 16 : AAG21_photospheric_element_zfrac(e_be) = 1.38d0 !± 0.09
1846 16 : AAG21_photospheric_element_zfrac(e_b ) = 2.70d0 !± 0.20
1847 16 : AAG21_photospheric_element_zfrac(e_c ) = 8.46d0 !± 0.04
1848 16 : AAG21_photospheric_element_zfrac(e_n ) = 7.83d0 !± 0.07
1849 16 : AAG21_photospheric_element_zfrac(e_o ) = 8.69d0 !± 0.04
1850 16 : AAG21_photospheric_element_zfrac(e_f ) = 4.40d0 !± 0.25
1851 16 : AAG21_photospheric_element_zfrac(e_ne) = 8.06d0 !± 0.05
1852 16 : AAG21_photospheric_element_zfrac(e_na) = 6.22d0 !± 0.03
1853 16 : AAG21_photospheric_element_zfrac(e_mg) = 7.55d0 !± 0.03
1854 16 : AAG21_photospheric_element_zfrac(e_al) = 6.43d0 !± 0.03
1855 16 : AAG21_photospheric_element_zfrac(e_si) = 7.51d0 !± 0.03
1856 16 : AAG21_photospheric_element_zfrac(e_p ) = 5.41d0 !± 0.03
1857 16 : AAG21_photospheric_element_zfrac(e_s ) = 7.12d0 !± 0.03
1858 16 : AAG21_photospheric_element_zfrac(e_cl) = 5.31d0 !± 0.20
1859 16 : AAG21_photospheric_element_zfrac(e_ar) = 6.38d0 !± 0.10
1860 16 : AAG21_photospheric_element_zfrac(e_k ) = 5.07d0 !± 0.03
1861 16 : AAG21_photospheric_element_zfrac(e_ca) = 6.30d0 !± 0.03
1862 16 : AAG21_photospheric_element_zfrac(e_sc) = 3.14d0 !± 0.04
1863 16 : AAG21_photospheric_element_zfrac(e_ti) = 4.97d0 !± 0.05
1864 16 : AAG21_photospheric_element_zfrac(e_v ) = 3.90d0 !± 0.08
1865 16 : AAG21_photospheric_element_zfrac(e_cr) = 5.62d0 !± 0.04
1866 16 : AAG21_photospheric_element_zfrac(e_mn) = 5.42d0 !± 0.06
1867 16 : AAG21_photospheric_element_zfrac(e_fe) = 7.46d0 !± 0.04
1868 16 : AAG21_photospheric_element_zfrac(e_co) = 4.94d0 !± 0.05
1869 16 : AAG21_photospheric_element_zfrac(e_ni) = 6.20d0 !± 0.04
1870 16 : AAG21_photospheric_element_zfrac(e_cu) = 4.18d0 !± 0.05
1871 16 : AAG21_photospheric_element_zfrac(e_zn) = 4.56d0 !± 0.05
1872 16 : AAG21_photospheric_element_zfrac(e_ga) = 3.02d0 !± 0.05
1873 16 : AAG21_photospheric_element_zfrac(e_ge) = 3.62d0 !± 0.10
1874 16 : AAG21_photospheric_element_zfrac(e_as) = 2.30d0 !± 0.04 meteorites
1875 16 : AAG21_photospheric_element_zfrac(e_se) = 3.34d0 !± 0.03 meteorites
1876 16 : AAG21_photospheric_element_zfrac(e_br) = 2.54d0 !± 0.06 meteorites
1877 16 : AAG21_photospheric_element_zfrac(e_kr) = 3.12d0 !± 0.10
1878 16 : AAG21_photospheric_element_zfrac(e_rb) = 2.32d0 !± 0.08
1879 16 : AAG21_photospheric_element_zfrac(e_sr) = 2.83d0 !± 0.06
1880 16 : AAG21_photospheric_element_zfrac(e_y ) = 2.21d0 !± 0.05
1881 16 : AAG21_photospheric_element_zfrac(e_zr) = 2.59d0 !± 0.04
1882 16 : AAG21_photospheric_element_zfrac(e_nb) = 1.47d0 !± 0.06
1883 16 : AAG21_photospheric_element_zfrac(e_mo) = 1.88d0 !± 0.09
1884 16 : AAG21_photospheric_element_zfrac(e_Ru) = 1.75d0 !± 0.08
1885 16 : AAG21_photospheric_element_zfrac(e_Rh) = 0.78d0 !± 0.11
1886 16 : AAG21_photospheric_element_zfrac(e_Pd) = 1.57d0 !± 0.10
1887 16 : AAG21_photospheric_element_zfrac(e_Ag) = 0.96d0 !± 0.10
1888 16 : AAG21_photospheric_element_zfrac(e_Cd) = 1.71d0 !± 0.03 meteorites
1889 16 : AAG21_photospheric_element_zfrac(e_In) = 0.80d0 !± 0.20
1890 16 : AAG21_photospheric_element_zfrac(e_Sn) = 2.02d0 !± 0.10
1891 16 : AAG21_photospheric_element_zfrac(e_Sb) = 1.01d0 !± 0.06 meteorites
1892 16 : AAG21_photospheric_element_zfrac(e_Te) = 2.18d0 !± 0.03 meteorites
1893 16 : AAG21_photospheric_element_zfrac(e_I ) = 1.55d0 !± 0.08 meteorites
1894 16 : AAG21_photospheric_element_zfrac(e_Xe) = 2.22d0 !± 0.05
1895 16 : AAG21_photospheric_element_zfrac(e_Cs) = 1.08d0 !± 0.03 meteorites
1896 16 : AAG21_photospheric_element_zfrac(e_Ba) = 2.27d0 !± 0.05
1897 16 : AAG21_photospheric_element_zfrac(e_La) = 1.11d0 !± 0.04
1898 16 : AAG21_photospheric_element_zfrac(e_Ce) = 1.58d0 !± 0.04
1899 16 : AAG21_photospheric_element_zfrac(e_Pr) = 0.75d0 !± 0.05
1900 16 : AAG21_photospheric_element_zfrac(e_Nd) = 1.42d0 !± 0.04
1901 16 : AAG21_photospheric_element_zfrac(e_Sm) = 0.95d0 !± 0.04
1902 16 : AAG21_photospheric_element_zfrac(e_Eu) = 0.52d0 !± 0.04
1903 16 : AAG21_photospheric_element_zfrac(e_Gd) = 1.08d0 !± 0.04
1904 16 : AAG21_photospheric_element_zfrac(e_Tb) = 0.31d0 !± 0.10
1905 16 : AAG21_photospheric_element_zfrac(e_Dy) = 1.10d0 !± 0.04
1906 16 : AAG21_photospheric_element_zfrac(e_Ho) = 0.48d0 !± 0.11
1907 16 : AAG21_photospheric_element_zfrac(e_Er) = 0.93d0 !± 0.05
1908 16 : AAG21_photospheric_element_zfrac(e_Tm) = 0.11d0 !± 0.04
1909 16 : AAG21_photospheric_element_zfrac(e_Yb) = 0.85d0 !± 0.11
1910 16 : AAG21_photospheric_element_zfrac(e_Lu) = 0.10d0 !± 0.09
1911 16 : AAG21_photospheric_element_zfrac(e_Hf) = 0.85d0 !± 0.05
1912 16 : AAG21_photospheric_element_zfrac(e_Ta) =-0.15d0 !± 0.04 meteorites
1913 16 : AAG21_photospheric_element_zfrac(e_W ) = 0.79d0 !± 0.11
1914 16 : AAG21_photospheric_element_zfrac(e_Re) = 0.26d0 !± 0.02 meteorites
1915 16 : AAG21_photospheric_element_zfrac(e_Os) = 1.35d0 !± 0.12
1916 16 : AAG21_photospheric_element_zfrac(e_Ir) = 1.32d0 !± 0.02 meteorites
1917 16 : AAG21_photospheric_element_zfrac(e_Pt) = 1.61d0 !± 0.02 meteorites
1918 16 : AAG21_photospheric_element_zfrac(e_Au) = 0.91d0 !± 0.12
1919 16 : AAG21_photospheric_element_zfrac(e_Hg) = 1.17d0 !± 0.18 meteorites
1920 16 : AAG21_photospheric_element_zfrac(e_Tl) = 0.92d0 !± 0.17
1921 16 : AAG21_photospheric_element_zfrac(e_Pb) = 1.95d0 !± 0.08
1922 16 : AAG21_photospheric_element_zfrac(e_Bi) = 0.65d0 !± 0.04 meteorites
1923 16 : AAG21_photospheric_element_zfrac(e_Th) = 0.03d0 !± 0.10
1924 16 : AAG21_photospheric_element_zfrac(e_U) =-0.54d0 !± 0.03 meteorites
1925 :
1926 : ! convert to fraction of Z by mass
1927 16 : z_sum = 0
1928 1456 : do i = e_li, e_u
1929 : AAG21_photospheric_element_zfrac(i) = &
1930 1440 : exp10(AAG21_photospheric_element_zfrac(i))*element_atomic_weight(i)
1931 1456 : z_sum = z_sum + AAG21_photospheric_element_zfrac(i)
1932 : end do
1933 1456 : do i = e_li, e_u
1934 1456 : AAG21_photospheric_element_zfrac(i) = AAG21_photospheric_element_zfrac(i) / z_sum
1935 : end do
1936 :
1937 16 : end subroutine init_AAG21_photospheric_data
1938 :
1939 16 : subroutine init_L09_data ! fraction by mass of total Z
1940 : ! Lodders 09
1941 : integer :: i
1942 16 : real(dp) :: z_sum
1943 : include 'formats'
1944 :
1945 16 : L09_element_zfrac(:) = 0
1946 :
1947 : ! mass fractions
1948 16 : L09_element_zfrac(e_li) = 1.054594933683d-08
1949 16 : L09_element_zfrac(e_be) = 1.5087555571d-10
1950 16 : L09_element_zfrac(e_b ) = 5.5633306761d-09
1951 16 : L09_element_zfrac(e_c ) = 0.002365544090848D0
1952 16 : L09_element_zfrac(e_n ) = 0.0008161934768925D0
1953 16 : L09_element_zfrac(e_o ) = 0.0068991682737478D0
1954 16 : L09_element_zfrac(e_f ) = 4.1844135602D-07
1955 16 : L09_element_zfrac(e_ne) = 0.0018162023028794D0
1956 16 : L09_element_zfrac(e_na) = 3.6352024324D-05
1957 16 : L09_element_zfrac(e_mg) = 0.000683514477408D0
1958 16 : L09_element_zfrac(e_al) = 6.2568980455D-05
1959 16 : L09_element_zfrac(e_si) = 0.000769722819625D0
1960 16 : L09_element_zfrac(e_p ) = 7.0479812061D-06
1961 16 : L09_element_zfrac(e_s ) = 0.000370123705833609D0
1962 16 : L09_element_zfrac(e_cl) = 5.0250763788D-06
1963 16 : L09_element_zfrac(e_ar) = 9.2220355099056D-05
1964 16 : L09_element_zfrac(e_k ) = 4.0297305059079D-06
1965 16 : L09_element_zfrac(e_ca) = 6.63135868398494D-05
1966 16 : L09_element_zfrac(e_sc) = 4.2402933957D-08
1967 16 : L09_element_zfrac(e_ti) = 3.24207135484D-06
1968 16 : L09_element_zfrac(e_v ) = 4.0049954611698D-07
1969 16 : L09_element_zfrac(e_cr) = 1.870484358513D-05
1970 16 : L09_element_zfrac(e_mn) = 1.3890521841D-05
1971 16 : L09_element_zfrac(e_fe) = 0.0012986862725768D0
1972 16 : L09_element_zfrac(e_co) = 3.7979113651D-06
1973 16 : L09_element_zfrac(e_ni) = 7.901833309378D-05
1974 16 : L09_element_zfrac(e_cu) = 9.4275308656D-07
1975 16 : L09_element_zfrac(e_zn) = 2.324135489336D-06
1976 16 : L09_element_zfrac(e_ga) = 6.9975797859D-08
1977 16 : L09_element_zfrac(e_ge) = 2.27918509226D-07
1978 16 : L09_element_zfrac(e_as) = 1.2531874861D-08
1979 16 : L09_element_zfrac(e_se) = 1.460613983153D-07
1980 16 : L09_element_zfrac(e_br) = 2.1909826069D-07
1981 16 : L09_element_zfrac(e_kr) = 1.2837104766341D-07
1982 16 : L09_element_zfrac(e_rb) = 1.69493949899D-08
1983 16 : L09_element_zfrac(e_sr) = 5.58119030984D-08
1984 16 : L09_element_zfrac(e_y ) = 1.1287452840D-08
1985 16 : L09_element_zfrac(e_zr) = 2.696974516888D-08
1986 16 : L09_element_zfrac(e_nb) = 1.9870212075D-09
1987 16 : L09_element_zfrac(e_mo) = 6.70675811282D-09
1988 16 : L09_element_zfrac(e_Ru) = 4.935148192679D-09
1989 16 : L09_element_zfrac(e_Rh) = 1.0439120240D-09
1990 16 : L09_element_zfrac(e_Pd) = 3.958834334296D-09
1991 16 : L09_element_zfrac(e_Ag) = 1.44909561511D-09
1992 16 : L09_element_zfrac(e_Cd) = 4.850725813727D-09
1993 16 : L09_element_zfrac(e_In) = 5.60277526584D-10
1994 16 : L09_element_zfrac(e_Sn) = 1.1748571050167D-08
1995 16 : L09_element_zfrac(e_Sb) = 1.04476117832D-09
1996 16 : L09_element_zfrac(e_Te) = 1.64381492798D-08
1997 16 : L09_element_zfrac(e_I ) = 3.8266730451D-09
1998 16 : L09_element_zfrac(e_Xe) = 1.9631791443838D-08
1999 16 : L09_element_zfrac(e_Cs) = 1.3516072159D-09
2000 16 : L09_element_zfrac(e_Ba) = 1.6848783894182D-08
2001 16 : L09_element_zfrac(e_La) = 1.7400268564D-09
2002 16 : L09_element_zfrac(e_Ce) = 4.5166246603333D-09
2003 16 : L09_element_zfrac(e_Pr) = 6.6431263199D-10
2004 16 : L09_element_zfrac(e_Nd) = 5.82366496835D-09
2005 16 : L09_element_zfrac(e_Sm) = 1.10017534847D-09
2006 16 : L09_element_zfrac(e_Eu) = 4.1023195079D-10
2007 16 : L09_element_zfrac(e_Gd) = 1.5505832574086D-09
2008 16 : L09_element_zfrac(e_Tb) = 2.7612856334D-10
2009 16 : L09_element_zfrac(e_Dy) = 1.79997715438722D-09
2010 16 : L09_element_zfrac(e_Ho) = 4.1129202414D-10
2011 16 : L09_element_zfrac(e_Er) = 1.2035968713737D-09
2012 16 : L09_element_zfrac(e_Tm) = 1.8794799164D-10
2013 16 : L09_element_zfrac(e_Yb) = 1.2153809425605D-09
2014 16 : L09_element_zfrac(e_Lu) = 1.826667993422D-10
2015 16 : L09_element_zfrac(e_Hf) = 7.619544268076D-10
2016 16 : L09_element_zfrac(e_Ta) = 1.04130101121661D-10
2017 16 : L09_element_zfrac(e_W ) = 6.9059258989518D-10
2018 16 : L09_element_zfrac(e_Re) = 2.9647265833D-10
2019 16 : L09_element_zfrac(e_Os) = 3.52828020217507D-09
2020 16 : L09_element_zfrac(e_Ir) = 3.5336600059D-09
2021 16 : L09_element_zfrac(e_Pt) = 6.8100262401586D-09
2022 16 : L09_element_zfrac(e_Au) = 1.0522666072D-09
2023 16 : L09_element_zfrac(e_Hg) = 2.5169209215851D-09
2024 16 : L09_element_zfrac(e_Tl) = 1.02465539439D-09
2025 16 : L09_element_zfrac(e_Pb) = 1.879940103275D-08
2026 16 : L09_element_zfrac(e_Bi) = 7.9004226175D-10
2027 16 : L09_element_zfrac(e_Th) = 2.7961831384D-10
2028 16 : L09_element_zfrac(e_U) = 1.546830543D-10
2029 :
2030 : ! convert from mass fraction to fraction of Z by mass
2031 1456 : z_sum = sum(L09_element_zfrac(e_li:e_u))
2032 1456 : do i = e_li, e_u
2033 1456 : L09_element_zfrac(i) = L09_element_zfrac(i) / z_sum
2034 : end do
2035 :
2036 16 : end subroutine init_L09_data
2037 :
2038 :
2039 16 : subroutine allocate_nuclide_data(d,n,ierr)
2040 : type(nuclide_data), intent(out) :: d
2041 : integer, intent(in) :: n
2042 : integer, intent(out) :: ierr
2043 16 : ierr = 0
2044 : allocate(d% name(n), d% W(n), d% Z(n), d% N(n), d% Z_plus_N(n), &
2045 : d% spin(N), d% binding_energy(n), d% Z53(n), &
2046 : d% isomeric_state(n), d% mass_excess(n), d% pfcn(npart,n), &
2047 16 : d% chem_id(n), d% nuclide(n), stat=ierr)
2048 16 : if (ierr /= 0) return
2049 16 : d% nnuclides = n
2050 : end subroutine allocate_nuclide_data
2051 :
2052 :
2053 1 : subroutine free_nuclide_data(n)
2054 : type(nuclide_data), intent(inout) :: n
2055 1 : if (associated(n% name)) &
2056 : deallocate( &
2057 0 : n% name, n% W, n% Z, n% N, n% Z_plus_N, n% spin, n% binding_energy, n% Z53, &
2058 1 : n% isomeric_state, n% mass_excess, n% pfcn, n% chem_id, n% nuclide)
2059 1 : n% nnuclides = 0
2060 1 : end subroutine free_nuclide_data
2061 :
2062 :
2063 1 : subroutine free_lodders03_table()
2064 : use utils_lib, only : integer_dict_free
2065 1 : deallocate(lodders03_tab6% isotopic_percent)
2066 1 : call integer_dict_free(lodders03_tab6% name_dict)
2067 1 : nullify(lodders03_tab6% name_dict)
2068 1 : end subroutine free_lodders03_table
2069 :
2070 :
2071 : ! returns the index of a particular nuclide in the full chem_isos set
2072 : ! returns nuclide_not_found if name not found
2073 98518 : function get_nuclide_index(nuclei) result(indx)
2074 : use utils_lib, only: integer_dict_lookup
2075 : character(len=*), intent(in) :: nuclei
2076 : integer :: indx, ierr
2077 49259 : if (.not. chem_has_been_initialized) then
2078 0 : write(*,*) 'must call chem_init before calling any other routine in chem_lib'
2079 : indx = nuclide_not_found
2080 0 : return
2081 : end if
2082 : ierr = 0
2083 49259 : call integer_dict_lookup(chem_isos_dict, nuclei, indx, ierr)
2084 49259 : if (ierr /= 0) indx = nuclide_not_found
2085 49259 : end function get_nuclide_index
2086 :
2087 :
2088 16 : subroutine set_some_isos
2089 16 : ih1 = get_nuclide_index('h1')
2090 16 : ih2 = get_nuclide_index('h2')
2091 16 : ih3 = get_nuclide_index('h3')
2092 16 : ihe3 = get_nuclide_index('he3')
2093 16 : ihe4 = get_nuclide_index('he4')
2094 16 : ili6 = get_nuclide_index('li6')
2095 16 : ili7 = get_nuclide_index('li7')
2096 16 : ili8 = get_nuclide_index('li8')
2097 16 : ibe7 = get_nuclide_index('be7')
2098 16 : ibe8 = get_nuclide_index('be8')
2099 16 : ibe9 = get_nuclide_index('be9')
2100 16 : ibe10 = get_nuclide_index('be10')
2101 16 : ibe11 = get_nuclide_index('be11')
2102 16 : ib8 = get_nuclide_index('b8')
2103 16 : ib10 = get_nuclide_index('b10')
2104 16 : ib11 = get_nuclide_index('b11')
2105 16 : ib12 = get_nuclide_index('b12')
2106 16 : ib13 = get_nuclide_index('b13')
2107 16 : ib14 = get_nuclide_index('b14')
2108 16 : ic9 = get_nuclide_index('c9')
2109 16 : ic10 = get_nuclide_index('c10')
2110 16 : ic11 = get_nuclide_index('c11')
2111 16 : ic12 = get_nuclide_index('c12')
2112 16 : ic13 = get_nuclide_index('c13')
2113 16 : ic14 = get_nuclide_index('c14')
2114 16 : ic15 = get_nuclide_index('c15')
2115 16 : ic16 = get_nuclide_index('c16')
2116 16 : in12 = get_nuclide_index('n12')
2117 16 : in13 = get_nuclide_index('n13')
2118 16 : in14 = get_nuclide_index('n14')
2119 16 : in15 = get_nuclide_index('n15')
2120 16 : in16 = get_nuclide_index('n16')
2121 16 : in17 = get_nuclide_index('n17')
2122 16 : in18 = get_nuclide_index('n18')
2123 16 : in19 = get_nuclide_index('n19')
2124 16 : in20 = get_nuclide_index('n20')
2125 16 : io13 = get_nuclide_index('o13')
2126 16 : io14 = get_nuclide_index('o14')
2127 16 : io15 = get_nuclide_index('o15')
2128 16 : io16 = get_nuclide_index('o16')
2129 16 : io17 = get_nuclide_index('o17')
2130 16 : io18 = get_nuclide_index('o18')
2131 16 : io19 = get_nuclide_index('o19')
2132 16 : io20 = get_nuclide_index('o20')
2133 16 : if15 = get_nuclide_index('f15')
2134 16 : if16 = get_nuclide_index('f16')
2135 16 : if17 = get_nuclide_index('f17')
2136 16 : if18 = get_nuclide_index('f18')
2137 16 : if19 = get_nuclide_index('f19')
2138 16 : if20 = get_nuclide_index('f20')
2139 16 : if21 = get_nuclide_index('f21')
2140 16 : if22 = get_nuclide_index('f22')
2141 16 : if23 = get_nuclide_index('f23')
2142 16 : if24 = get_nuclide_index('f24')
2143 16 : ine17 = get_nuclide_index('ne17')
2144 16 : ine18 = get_nuclide_index('ne18')
2145 16 : ine19 = get_nuclide_index('ne19')
2146 16 : ine20 = get_nuclide_index('ne20')
2147 16 : ine21 = get_nuclide_index('ne21')
2148 16 : ine22 = get_nuclide_index('ne22')
2149 16 : ine23 = get_nuclide_index('ne23')
2150 16 : ine24 = get_nuclide_index('ne24')
2151 16 : ine25 = get_nuclide_index('ne25')
2152 16 : ine26 = get_nuclide_index('ne26')
2153 16 : ine27 = get_nuclide_index('ne27')
2154 16 : ine28 = get_nuclide_index('ne28')
2155 16 : ina20 = get_nuclide_index('na20')
2156 16 : ina21 = get_nuclide_index('na21')
2157 16 : ina22 = get_nuclide_index('na22')
2158 16 : ina23 = get_nuclide_index('na23')
2159 16 : ina24 = get_nuclide_index('na24')
2160 16 : ina25 = get_nuclide_index('na25')
2161 16 : ina26 = get_nuclide_index('na26')
2162 16 : ina27 = get_nuclide_index('na27')
2163 16 : ina28 = get_nuclide_index('na28')
2164 16 : ina29 = get_nuclide_index('na29')
2165 16 : ina30 = get_nuclide_index('na30')
2166 16 : ina31 = get_nuclide_index('na31')
2167 16 : img20 = get_nuclide_index('mg20')
2168 16 : img21 = get_nuclide_index('mg21')
2169 16 : img22 = get_nuclide_index('mg22')
2170 16 : img23 = get_nuclide_index('mg23')
2171 16 : img24 = get_nuclide_index('mg24')
2172 16 : img25 = get_nuclide_index('mg25')
2173 16 : img26 = get_nuclide_index('mg26')
2174 16 : img27 = get_nuclide_index('mg27')
2175 16 : img28 = get_nuclide_index('mg28')
2176 16 : img29 = get_nuclide_index('mg29')
2177 16 : img30 = get_nuclide_index('mg30')
2178 16 : img31 = get_nuclide_index('mg31')
2179 16 : img32 = get_nuclide_index('mg32')
2180 16 : img33 = get_nuclide_index('mg33')
2181 16 : ial22 = get_nuclide_index('al22')
2182 16 : ial23 = get_nuclide_index('al23')
2183 16 : ial24 = get_nuclide_index('al24')
2184 16 : ial25 = get_nuclide_index('al25')
2185 16 : ial26 = get_nuclide_index('al26')
2186 16 : ial27 = get_nuclide_index('al27')
2187 16 : ial28 = get_nuclide_index('al28')
2188 16 : ial29 = get_nuclide_index('al29')
2189 16 : ial30 = get_nuclide_index('al30')
2190 16 : ial31 = get_nuclide_index('al31')
2191 16 : ial32 = get_nuclide_index('al32')
2192 16 : ial33 = get_nuclide_index('al33')
2193 16 : ial34 = get_nuclide_index('al34')
2194 16 : ial35 = get_nuclide_index('al35')
2195 16 : isi22 = get_nuclide_index('si22')
2196 16 : isi23 = get_nuclide_index('si23')
2197 16 : isi24 = get_nuclide_index('si24')
2198 16 : isi25 = get_nuclide_index('si25')
2199 16 : isi26 = get_nuclide_index('si26')
2200 16 : isi27 = get_nuclide_index('si27')
2201 16 : isi28 = get_nuclide_index('si28')
2202 16 : isi29 = get_nuclide_index('si29')
2203 16 : isi30 = get_nuclide_index('si30')
2204 16 : isi31 = get_nuclide_index('si31')
2205 16 : isi32 = get_nuclide_index('si32')
2206 16 : isi33 = get_nuclide_index('si33')
2207 16 : isi34 = get_nuclide_index('si34')
2208 16 : isi35 = get_nuclide_index('si35')
2209 16 : isi36 = get_nuclide_index('si36')
2210 16 : isi37 = get_nuclide_index('si37')
2211 16 : isi38 = get_nuclide_index('si38')
2212 16 : ip26 = get_nuclide_index('p26')
2213 16 : ip27 = get_nuclide_index('p27')
2214 16 : ip28 = get_nuclide_index('p28')
2215 16 : ip29 = get_nuclide_index('p29')
2216 16 : ip30 = get_nuclide_index('p30')
2217 16 : ip31 = get_nuclide_index('p31')
2218 16 : ip32 = get_nuclide_index('p32')
2219 16 : ip33 = get_nuclide_index('p33')
2220 16 : ip34 = get_nuclide_index('p34')
2221 16 : ip35 = get_nuclide_index('p35')
2222 16 : ip36 = get_nuclide_index('p36')
2223 16 : ip37 = get_nuclide_index('p37')
2224 16 : ip38 = get_nuclide_index('p38')
2225 16 : ip39 = get_nuclide_index('p39')
2226 16 : ip40 = get_nuclide_index('p40')
2227 16 : is27 = get_nuclide_index('s27')
2228 16 : is28 = get_nuclide_index('s28')
2229 16 : is29 = get_nuclide_index('s29')
2230 16 : is30 = get_nuclide_index('s30')
2231 16 : is31 = get_nuclide_index('s31')
2232 16 : is32 = get_nuclide_index('s32')
2233 16 : is33 = get_nuclide_index('s33')
2234 16 : is34 = get_nuclide_index('s34')
2235 16 : is35 = get_nuclide_index('s35')
2236 16 : is36 = get_nuclide_index('s36')
2237 16 : is37 = get_nuclide_index('s37')
2238 16 : is38 = get_nuclide_index('s38')
2239 16 : is39 = get_nuclide_index('s39')
2240 16 : is40 = get_nuclide_index('s40')
2241 16 : is41 = get_nuclide_index('s41')
2242 16 : is42 = get_nuclide_index('s42')
2243 16 : icl31 = get_nuclide_index('cl31')
2244 16 : icl32 = get_nuclide_index('cl32')
2245 16 : icl33 = get_nuclide_index('cl33')
2246 16 : icl34 = get_nuclide_index('cl34')
2247 16 : icl35 = get_nuclide_index('cl35')
2248 16 : icl36 = get_nuclide_index('cl36')
2249 16 : icl37 = get_nuclide_index('cl37')
2250 16 : icl38 = get_nuclide_index('cl38')
2251 16 : icl39 = get_nuclide_index('cl39')
2252 16 : icl40 = get_nuclide_index('cl40')
2253 16 : icl41 = get_nuclide_index('cl41')
2254 16 : icl42 = get_nuclide_index('cl42')
2255 16 : icl43 = get_nuclide_index('cl43')
2256 16 : icl44 = get_nuclide_index('cl44')
2257 16 : iar31 = get_nuclide_index('ar31')
2258 16 : iar32 = get_nuclide_index('ar32')
2259 16 : iar33 = get_nuclide_index('ar33')
2260 16 : iar34 = get_nuclide_index('ar34')
2261 16 : iar35 = get_nuclide_index('ar35')
2262 16 : iar36 = get_nuclide_index('ar36')
2263 16 : iar37 = get_nuclide_index('ar37')
2264 16 : iar38 = get_nuclide_index('ar38')
2265 16 : iar39 = get_nuclide_index('ar39')
2266 16 : iar40 = get_nuclide_index('ar40')
2267 16 : iar41 = get_nuclide_index('ar41')
2268 16 : iar42 = get_nuclide_index('ar42')
2269 16 : iar43 = get_nuclide_index('ar43')
2270 16 : iar44 = get_nuclide_index('ar44')
2271 16 : iar45 = get_nuclide_index('ar45')
2272 16 : iar46 = get_nuclide_index('ar46')
2273 16 : iar47 = get_nuclide_index('ar47')
2274 16 : ik35 = get_nuclide_index('k35')
2275 16 : ik36 = get_nuclide_index('k36')
2276 16 : ik37 = get_nuclide_index('k37')
2277 16 : ik38 = get_nuclide_index('k38')
2278 16 : ik39 = get_nuclide_index('k39')
2279 16 : ik40 = get_nuclide_index('k40')
2280 16 : ik41 = get_nuclide_index('k41')
2281 16 : ik42 = get_nuclide_index('k42')
2282 16 : ik43 = get_nuclide_index('k43')
2283 16 : ik44 = get_nuclide_index('k44')
2284 16 : ik45 = get_nuclide_index('k45')
2285 16 : ik46 = get_nuclide_index('k46')
2286 16 : ik47 = get_nuclide_index('k47')
2287 16 : ica35 = get_nuclide_index('ca35')
2288 16 : ica36 = get_nuclide_index('ca36')
2289 16 : ica37 = get_nuclide_index('ca37')
2290 16 : ica38 = get_nuclide_index('ca38')
2291 16 : ica39 = get_nuclide_index('ca39')
2292 16 : ica40 = get_nuclide_index('ca40')
2293 16 : ica41 = get_nuclide_index('ca41')
2294 16 : ica42 = get_nuclide_index('ca42')
2295 16 : ica43 = get_nuclide_index('ca43')
2296 16 : ica44 = get_nuclide_index('ca44')
2297 16 : ica45 = get_nuclide_index('ca45')
2298 16 : ica46 = get_nuclide_index('ca46')
2299 16 : ica47 = get_nuclide_index('ca47')
2300 16 : ica48 = get_nuclide_index('ca48')
2301 16 : ica49 = get_nuclide_index('ca49')
2302 16 : ica50 = get_nuclide_index('ca50')
2303 16 : ica51 = get_nuclide_index('ca51')
2304 16 : ica52 = get_nuclide_index('ca52')
2305 16 : ica53 = get_nuclide_index('ca53')
2306 16 : isc40 = get_nuclide_index('sc40')
2307 16 : isc41 = get_nuclide_index('sc41')
2308 16 : isc42 = get_nuclide_index('sc42')
2309 16 : isc43 = get_nuclide_index('sc43')
2310 16 : isc44 = get_nuclide_index('sc44')
2311 16 : isc45 = get_nuclide_index('sc45')
2312 16 : isc46 = get_nuclide_index('sc46')
2313 16 : isc47 = get_nuclide_index('sc47')
2314 16 : isc48 = get_nuclide_index('sc48')
2315 16 : isc49 = get_nuclide_index('sc49')
2316 16 : isc50 = get_nuclide_index('sc50')
2317 16 : isc51 = get_nuclide_index('sc51')
2318 16 : isc52 = get_nuclide_index('sc52')
2319 16 : isc53 = get_nuclide_index('sc53')
2320 16 : iti39 = get_nuclide_index('ti39')
2321 16 : iti40 = get_nuclide_index('ti40')
2322 16 : iti41 = get_nuclide_index('ti41')
2323 16 : iti42 = get_nuclide_index('ti42')
2324 16 : iti43 = get_nuclide_index('ti43')
2325 16 : iti44 = get_nuclide_index('ti44')
2326 16 : iti45 = get_nuclide_index('ti45')
2327 16 : iti46 = get_nuclide_index('ti46')
2328 16 : iti47 = get_nuclide_index('ti47')
2329 16 : iti48 = get_nuclide_index('ti48')
2330 16 : iti49 = get_nuclide_index('ti49')
2331 16 : iti50 = get_nuclide_index('ti50')
2332 16 : iti51 = get_nuclide_index('ti51')
2333 16 : iti52 = get_nuclide_index('ti52')
2334 16 : iti53 = get_nuclide_index('ti53')
2335 16 : iti54 = get_nuclide_index('ti54')
2336 16 : iti55 = get_nuclide_index('ti55')
2337 16 : iv43 = get_nuclide_index('v43')
2338 16 : iv44 = get_nuclide_index('v44')
2339 16 : iv45 = get_nuclide_index('v45')
2340 16 : iv46 = get_nuclide_index('v46')
2341 16 : iv47 = get_nuclide_index('v47')
2342 16 : iv48 = get_nuclide_index('v48')
2343 16 : iv49 = get_nuclide_index('v49')
2344 16 : iv50 = get_nuclide_index('v50')
2345 16 : iv51 = get_nuclide_index('v51')
2346 16 : iv52 = get_nuclide_index('v52')
2347 16 : iv53 = get_nuclide_index('v53')
2348 16 : iv54 = get_nuclide_index('v54')
2349 16 : iv55 = get_nuclide_index('v55')
2350 16 : iv56 = get_nuclide_index('v56')
2351 16 : iv57 = get_nuclide_index('v57')
2352 16 : icr43 = get_nuclide_index('cr43')
2353 16 : icr44 = get_nuclide_index('cr44')
2354 16 : icr45 = get_nuclide_index('cr45')
2355 16 : icr46 = get_nuclide_index('cr46')
2356 16 : icr47 = get_nuclide_index('cr47')
2357 16 : icr48 = get_nuclide_index('cr48')
2358 16 : icr49 = get_nuclide_index('cr49')
2359 16 : icr50 = get_nuclide_index('cr50')
2360 16 : icr51 = get_nuclide_index('cr51')
2361 16 : icr52 = get_nuclide_index('cr52')
2362 16 : icr53 = get_nuclide_index('cr53')
2363 16 : icr54 = get_nuclide_index('cr54')
2364 16 : icr55 = get_nuclide_index('cr55')
2365 16 : icr56 = get_nuclide_index('cr56')
2366 16 : icr57 = get_nuclide_index('cr57')
2367 16 : icr58 = get_nuclide_index('cr58')
2368 16 : icr59 = get_nuclide_index('cr59')
2369 16 : icr60 = get_nuclide_index('cr60')
2370 16 : icr61 = get_nuclide_index('cr61')
2371 16 : icr62 = get_nuclide_index('cr62')
2372 16 : icr63 = get_nuclide_index('cr63')
2373 16 : icr64 = get_nuclide_index('cr64')
2374 16 : icr65 = get_nuclide_index('cr65')
2375 16 : icr66 = get_nuclide_index('cr66')
2376 16 : imn46 = get_nuclide_index('mn46')
2377 16 : imn47 = get_nuclide_index('mn47')
2378 16 : imn48 = get_nuclide_index('mn48')
2379 16 : imn49 = get_nuclide_index('mn49')
2380 16 : imn50 = get_nuclide_index('mn50')
2381 16 : imn51 = get_nuclide_index('mn51')
2382 16 : imn52 = get_nuclide_index('mn52')
2383 16 : imn53 = get_nuclide_index('mn53')
2384 16 : imn54 = get_nuclide_index('mn54')
2385 16 : imn55 = get_nuclide_index('mn55')
2386 16 : imn56 = get_nuclide_index('mn56')
2387 16 : imn57 = get_nuclide_index('mn57')
2388 16 : imn58 = get_nuclide_index('mn58')
2389 16 : imn59 = get_nuclide_index('mn59')
2390 16 : imn60 = get_nuclide_index('mn60')
2391 16 : imn61 = get_nuclide_index('mn61')
2392 16 : imn62 = get_nuclide_index('mn62')
2393 16 : imn63 = get_nuclide_index('mn63')
2394 16 : ife46 = get_nuclide_index('fe46')
2395 16 : ife47 = get_nuclide_index('fe47')
2396 16 : ife48 = get_nuclide_index('fe48')
2397 16 : ife49 = get_nuclide_index('fe49')
2398 16 : ife50 = get_nuclide_index('fe50')
2399 16 : ife51 = get_nuclide_index('fe51')
2400 16 : ife52 = get_nuclide_index('fe52')
2401 16 : ife53 = get_nuclide_index('fe53')
2402 16 : ife54 = get_nuclide_index('fe54')
2403 16 : ife55 = get_nuclide_index('fe55')
2404 16 : ife56 = get_nuclide_index('fe56')
2405 16 : ife57 = get_nuclide_index('fe57')
2406 16 : ife58 = get_nuclide_index('fe58')
2407 16 : ife59 = get_nuclide_index('fe59')
2408 16 : ife60 = get_nuclide_index('fe60')
2409 16 : ife61 = get_nuclide_index('fe61')
2410 16 : ife62 = get_nuclide_index('fe62')
2411 16 : ife63 = get_nuclide_index('fe63')
2412 16 : ife64 = get_nuclide_index('fe64')
2413 16 : ife65 = get_nuclide_index('fe65')
2414 16 : ife66 = get_nuclide_index('fe66')
2415 16 : ife68 = get_nuclide_index('fe68')
2416 16 : ico50 = get_nuclide_index('co50')
2417 16 : ico51 = get_nuclide_index('co51')
2418 16 : ico52 = get_nuclide_index('co52')
2419 16 : ico53 = get_nuclide_index('co53')
2420 16 : ico54 = get_nuclide_index('co54')
2421 16 : ico55 = get_nuclide_index('co55')
2422 16 : ico56 = get_nuclide_index('co56')
2423 16 : ico57 = get_nuclide_index('co57')
2424 16 : ico58 = get_nuclide_index('co58')
2425 16 : ico59 = get_nuclide_index('co59')
2426 16 : ico60 = get_nuclide_index('co60')
2427 16 : ico61 = get_nuclide_index('co61')
2428 16 : ico62 = get_nuclide_index('co62')
2429 16 : ico63 = get_nuclide_index('co63')
2430 16 : ico64 = get_nuclide_index('co64')
2431 16 : ico65 = get_nuclide_index('co65')
2432 16 : ico66 = get_nuclide_index('co66')
2433 16 : ico67 = get_nuclide_index('co67')
2434 16 : ini50 = get_nuclide_index('ni50')
2435 16 : ini51 = get_nuclide_index('ni51')
2436 16 : ini52 = get_nuclide_index('ni52')
2437 16 : ini53 = get_nuclide_index('ni53')
2438 16 : ini54 = get_nuclide_index('ni54')
2439 16 : ini55 = get_nuclide_index('ni55')
2440 16 : ini56 = get_nuclide_index('ni56')
2441 16 : ini57 = get_nuclide_index('ni57')
2442 16 : ini58 = get_nuclide_index('ni58')
2443 16 : ini59 = get_nuclide_index('ni59')
2444 16 : ini60 = get_nuclide_index('ni60')
2445 16 : ini61 = get_nuclide_index('ni61')
2446 16 : ini62 = get_nuclide_index('ni62')
2447 16 : ini63 = get_nuclide_index('ni63')
2448 16 : ini64 = get_nuclide_index('ni64')
2449 16 : ini65 = get_nuclide_index('ni65')
2450 16 : ini66 = get_nuclide_index('ni66')
2451 16 : ini67 = get_nuclide_index('ni67')
2452 16 : ini68 = get_nuclide_index('ni68')
2453 16 : ini69 = get_nuclide_index('ni69')
2454 16 : ini70 = get_nuclide_index('ni70')
2455 16 : ini71 = get_nuclide_index('ni71')
2456 16 : ini72 = get_nuclide_index('ni72')
2457 16 : ini73 = get_nuclide_index('ni73')
2458 16 : icu56 = get_nuclide_index('cu56')
2459 16 : icu57 = get_nuclide_index('cu57')
2460 16 : icu58 = get_nuclide_index('cu58')
2461 16 : icu59 = get_nuclide_index('cu59')
2462 16 : icu60 = get_nuclide_index('cu60')
2463 16 : icu61 = get_nuclide_index('cu61')
2464 16 : icu62 = get_nuclide_index('cu62')
2465 16 : icu63 = get_nuclide_index('cu63')
2466 16 : icu64 = get_nuclide_index('cu64')
2467 16 : icu65 = get_nuclide_index('cu65')
2468 16 : icu66 = get_nuclide_index('cu66')
2469 16 : icu67 = get_nuclide_index('cu67')
2470 16 : icu68 = get_nuclide_index('cu68')
2471 16 : icu69 = get_nuclide_index('cu69')
2472 16 : icu70 = get_nuclide_index('cu70')
2473 16 : icu71 = get_nuclide_index('cu71')
2474 16 : icu72 = get_nuclide_index('cu72')
2475 16 : izn55 = get_nuclide_index('zn55')
2476 16 : izn56 = get_nuclide_index('zn56')
2477 16 : izn57 = get_nuclide_index('zn57')
2478 16 : izn58 = get_nuclide_index('zn58')
2479 16 : izn59 = get_nuclide_index('zn59')
2480 16 : izn60 = get_nuclide_index('zn60')
2481 16 : izn61 = get_nuclide_index('zn61')
2482 16 : izn62 = get_nuclide_index('zn62')
2483 16 : izn63 = get_nuclide_index('zn63')
2484 16 : izn64 = get_nuclide_index('zn64')
2485 16 : izn65 = get_nuclide_index('zn65')
2486 16 : izn66 = get_nuclide_index('zn66')
2487 16 : izn67 = get_nuclide_index('zn67')
2488 16 : izn68 = get_nuclide_index('zn68')
2489 16 : izn69 = get_nuclide_index('zn69')
2490 16 : izn70 = get_nuclide_index('zn70')
2491 16 : izn71 = get_nuclide_index('zn71')
2492 16 : izn72 = get_nuclide_index('zn72')
2493 16 : izn73 = get_nuclide_index('zn73')
2494 16 : izn74 = get_nuclide_index('zn74')
2495 :
2496 16 : iga60 = get_nuclide_index('ga60')
2497 16 : iga61 = get_nuclide_index('ga61')
2498 16 : iga62 = get_nuclide_index('ga62')
2499 16 : iga63 = get_nuclide_index('ga63')
2500 16 : iga64 = get_nuclide_index('ga64')
2501 16 : iga65 = get_nuclide_index('ga65')
2502 16 : iga66 = get_nuclide_index('ga66')
2503 16 : iga67 = get_nuclide_index('ga67')
2504 16 : iga68 = get_nuclide_index('ga68')
2505 16 : iga69 = get_nuclide_index('ga69')
2506 16 : iga70 = get_nuclide_index('ga70')
2507 16 : iga71 = get_nuclide_index('ga71')
2508 16 : iga72 = get_nuclide_index('ga72')
2509 16 : iga73 = get_nuclide_index('ga73')
2510 16 : iga74 = get_nuclide_index('ga74')
2511 16 : iga75 = get_nuclide_index('ga75')
2512 :
2513 16 : ige59 = get_nuclide_index('ge59')
2514 16 : ige60 = get_nuclide_index('ge60')
2515 16 : ige61 = get_nuclide_index('ge61')
2516 16 : ige62 = get_nuclide_index('ge62')
2517 16 : ige63 = get_nuclide_index('ge63')
2518 16 : ige64 = get_nuclide_index('ge64')
2519 16 : ige65 = get_nuclide_index('ge65')
2520 16 : ige66 = get_nuclide_index('ge66')
2521 16 : ige67 = get_nuclide_index('ge67')
2522 16 : ige68 = get_nuclide_index('ge68')
2523 16 : ige69 = get_nuclide_index('ge69')
2524 16 : ige70 = get_nuclide_index('ge70')
2525 16 : ige71 = get_nuclide_index('ge71')
2526 16 : ige72 = get_nuclide_index('ge72')
2527 16 : ige73 = get_nuclide_index('ge73')
2528 16 : ige74 = get_nuclide_index('ge74')
2529 16 : ige75 = get_nuclide_index('ge75')
2530 16 : ige76 = get_nuclide_index('ge76')
2531 :
2532 16 : ias71 = get_nuclide_index('as71')
2533 16 : ias72 = get_nuclide_index('as72')
2534 16 : ias73 = get_nuclide_index('as73')
2535 16 : ias74 = get_nuclide_index('as74')
2536 16 : ias75 = get_nuclide_index('as75')
2537 16 : ias76 = get_nuclide_index('as76')
2538 16 : ias77 = get_nuclide_index('as77')
2539 16 : ias78 = get_nuclide_index('as78')
2540 16 : ias79 = get_nuclide_index('as79')
2541 :
2542 16 : ise68 = get_nuclide_index('se68')
2543 16 : ise69 = get_nuclide_index('se69')
2544 16 : ise70 = get_nuclide_index('se70')
2545 16 : ise71 = get_nuclide_index('se71')
2546 16 : ise72 = get_nuclide_index('se72')
2547 16 : ise73 = get_nuclide_index('se73')
2548 16 : ise74 = get_nuclide_index('se74')
2549 16 : ise75 = get_nuclide_index('se75')
2550 16 : ise76 = get_nuclide_index('se76')
2551 :
2552 16 : ikr70 = get_nuclide_index('kr70')
2553 16 : ikr72 = get_nuclide_index('kr72')
2554 16 : isr74 = get_nuclide_index('sr74')
2555 16 : isr76 = get_nuclide_index('sr76')
2556 16 : izr77 = get_nuclide_index('zr77')
2557 16 : izr80 = get_nuclide_index('zr80')
2558 16 : imo82 = get_nuclide_index('mo82')
2559 16 : imo84 = get_nuclide_index('mo84')
2560 16 : iru88 = get_nuclide_index('ru88')
2561 16 : ipd92 = get_nuclide_index('pd92')
2562 16 : isr75 = get_nuclide_index('sr75')
2563 16 : iru86 = get_nuclide_index('ru86')
2564 16 : iru87 = get_nuclide_index('ru87')
2565 16 : ipd89 = get_nuclide_index('pd89')
2566 16 : ipd91 = get_nuclide_index('pd91')
2567 16 : icd93 = get_nuclide_index('cd93')
2568 16 : icd96 = get_nuclide_index('cd96')
2569 16 : isn98 = get_nuclide_index('sn98')
2570 16 : isn100 = get_nuclide_index('sn100')
2571 16 : isn102 = get_nuclide_index('sn102')
2572 16 : isn104 = get_nuclide_index('sn104')
2573 16 : ineut = get_nuclide_index('neut')
2574 16 : iprot = get_nuclide_index('prot')
2575 :
2576 16 : end subroutine set_some_isos
2577 :
2578 :
2579 1590 : integer function category_id(cname)
2580 : character (len=*), intent(in) :: cname
2581 : ! returns id for the category if there is a matching name
2582 : ! returns 0 otherwise.
2583 : integer :: i, len
2584 : character (len=maxlen_category_name) :: nam
2585 1590 : len = len_trim(cname)
2586 27030 : do i = 1, maxlen_category_name
2587 27030 : if (i <= len) then
2588 8850 : nam(i:i) = cname(i:i)
2589 : else
2590 16590 : nam(i:i) = ' '
2591 : end if
2592 : end do
2593 23275 : do i = 1, num_categories
2594 23275 : if (category_name(i)==nam) then
2595 1590 : category_id = i
2596 1590 : return
2597 : end if
2598 : end do
2599 1590 : category_id = 0
2600 : end function category_id
2601 :
2602 :
2603 16 : subroutine set_category_names
2604 : integer :: i
2605 400 : category_name(:) = ''
2606 :
2607 16 : category_name(ipp) = 'pp'
2608 16 : category_name(icno) = 'cno'
2609 16 : category_name(i3alf) = 'tri_alpha'
2610 :
2611 16 : category_name(i_burn_c) = 'c_alpha'
2612 16 : category_name(i_burn_n) = 'n_alpha'
2613 16 : category_name(i_burn_o) = 'o_alpha'
2614 16 : category_name(i_burn_ne) = 'ne_alpha'
2615 16 : category_name(i_burn_na) = 'na_alpha'
2616 16 : category_name(i_burn_mg) = 'mg_alpha'
2617 16 : category_name(i_burn_si) = 'si_alpha'
2618 16 : category_name(i_burn_s) = 's_alpha'
2619 16 : category_name(i_burn_ar) = 'ar_alpha'
2620 16 : category_name(i_burn_ca) = 'ca_alpha'
2621 16 : category_name(i_burn_ti) = 'ti_alpha'
2622 16 : category_name(i_burn_cr) = 'cr_alpha'
2623 16 : category_name(i_burn_fe) = 'fe_co_ni'
2624 :
2625 16 : category_name(icc) = 'c12_c12'
2626 16 : category_name(ico) = 'c12_o16'
2627 16 : category_name(ioo) = 'o16_o16'
2628 :
2629 16 : category_name(iphoto) = 'photo'
2630 16 : category_name(ipnhe4) = 'pnhe4'
2631 :
2632 16 : category_name(i_ni56_co56) = 'ni56_co56'
2633 16 : category_name(i_co56_fe56) = 'co56_fe56'
2634 :
2635 16 : category_name(iother) = 'other'
2636 :
2637 16 : i=1 ! write it this way to avoid stupid compiler warning.
2638 16 : if (len_trim(category_name(i)) == 0) then
2639 0 : write(*,*) 'missing name for category', i
2640 0 : flush(6)
2641 0 : call mesa_error(__FILE__,__LINE__,'set_category_names')
2642 : end if
2643 384 : do i=2,num_categories
2644 384 : if (len_trim(category_name(i)) == 0) then
2645 0 : write(*,*) 'missing name for category', i
2646 0 : write(*,*) 'following ' // trim(category_name(i-1))
2647 0 : flush(6)
2648 0 : call mesa_error(__FILE__,__LINE__,'set_category_names')
2649 : end if
2650 : end do
2651 :
2652 16 : end subroutine set_category_names
2653 :
2654 0 : end module chem_def
|