9
9
from scipy .interpolate import CubicSpline , interp1d , PchipInterpolator
10
10
11
11
from fiasco import proton_electron_ratio
12
- from fiasco .base import ContinuumBase , IonBase
12
+ from fiasco .base import IonBase
13
13
from fiasco .collections import IonCollection
14
14
from fiasco .gaunt import GauntFactor
15
15
from fiasco .levels import Level , Transitions
24
24
__all__ = ['Ion' ]
25
25
26
26
27
- class Ion (IonBase , ContinuumBase ):
27
+ class Ion (IonBase ):
28
28
"""
29
29
Class for representing a CHIANTI ion.
30
30
@@ -52,11 +52,14 @@ class Ion(IonBase, ContinuumBase):
52
52
"""
53
53
54
54
@u .quantity_input
55
- def __init__ (self , ion_name , temperature : u .K ,
56
- abundance = 'sun_coronal_1992_feldman_ext' ,
57
- ionization_fraction = 'chianti' ,
58
- ionization_potential = 'chianti' ,
59
- * args , ** kwargs ):
55
+ def __init__ (self ,
56
+ ion_name ,
57
+ temperature : u .K ,
58
+ abundance = 'sun_coronal_1992_feldman_ext' ,
59
+ ionization_fraction = 'chianti' ,
60
+ ionization_potential = 'chianti' ,
61
+ * args ,
62
+ ** kwargs ):
60
63
super ().__init__ (ion_name , * args , ** kwargs )
61
64
self .temperature = np .atleast_1d (temperature )
62
65
self ._dset_names = {}
@@ -163,6 +166,14 @@ def _has_dataset(self, dset_name):
163
166
else :
164
167
return True
165
168
169
+ @property
170
+ @u .quantity_input
171
+ def thermal_energy (self ) -> u .erg :
172
+ """
173
+ Thermal energy, :math:`k_BT`, as a function of temperature.
174
+ """
175
+ return self .temperature .to ('erg' , equivalencies = u .equivalencies .temperature_energy ())
176
+
166
177
def next_ion (self ):
167
178
"""
168
179
Return an `~fiasco.Ion` instance with the next highest ionization stage.
@@ -373,7 +384,7 @@ def effective_collision_strength(self) -> u.dimensionless_unscaled:
373
384
--------
374
385
fiasco.util.burgess_tully_descale : Descale and interpolate :math:`\Upsilon`.
375
386
"""
376
- kBTE = np .outer (const . k_B * self .temperature , 1.0 / self ._scups ['delta_energy' ])
387
+ kBTE = np .outer (self .thermal_energy , 1.0 / self ._scups ['delta_energy' ])
377
388
upsilon = burgess_tully_descale (self ._scups ['bt_t' ],
378
389
self ._scups ['bt_upsilon' ],
379
390
kBTE .T ,
@@ -406,10 +417,10 @@ def electron_collision_deexcitation_rate(self) -> u.cm**3 / u.s:
406
417
electron_collision_excitation_rate : Excitation rate due to collisions
407
418
effective_collision_strength : Maxwellian-averaged collision strength, :math:`\Upsilon`
408
419
"""
409
- c = ( const .h ** 2 ) / (( 2. * np .pi * const .m_e )** (1.5 ) * np . sqrt ( const . k_B ) )
420
+ c = const .h ** 2 / (2. * np .pi * const .m_e )** (1.5 )
410
421
upsilon = self .effective_collision_strength
411
422
omega_upper = 2. * self ._elvlc ['J' ][self ._scups ['upper_level' ] - 1 ] + 1.
412
- return c * upsilon / np .sqrt (self .temperature [:, np .newaxis ]) / omega_upper
423
+ return c * upsilon / np .sqrt (self .thermal_energy [:, np .newaxis ]) / omega_upper
413
424
414
425
@cached_property
415
426
@needs_dataset ('elvlc' , 'scups' )
@@ -439,7 +450,7 @@ def electron_collision_excitation_rate(self) -> u.cm**3 / u.s:
439
450
"""
440
451
omega_upper = 2. * self ._elvlc ['J' ][self ._scups ['upper_level' ] - 1 ] + 1.
441
452
omega_lower = 2. * self ._elvlc ['J' ][self ._scups ['lower_level' ] - 1 ] + 1.
442
- kBTE = np .outer (1. / const . k_B / self .temperature , self ._scups ['delta_energy' ])
453
+ kBTE = np .outer (1. / self .thermal_energy , self ._scups ['delta_energy' ])
443
454
return omega_upper / omega_lower * self .electron_collision_deexcitation_rate * np .exp (- kBTE )
444
455
445
456
@cached_property
@@ -461,7 +472,7 @@ def proton_collision_excitation_rate(self) -> u.cm**3 / u.s:
461
472
# Create scaled temperature--these are not stored in the file
462
473
bt_t = [np .linspace (0 , 1 , ups .shape [0 ]) for ups in self ._psplups ['bt_rate' ]]
463
474
# Get excitation rates directly from scaled data
464
- kBTE = np .outer (const . k_B * self .temperature , 1.0 / self ._psplups ['delta_energy' ])
475
+ kBTE = np .outer (self .thermal_energy , 1.0 / self ._psplups ['delta_energy' ])
465
476
ex_rate = burgess_tully_descale (bt_t ,
466
477
self ._psplups ['bt_rate' ],
467
478
kBTE .T ,
@@ -494,7 +505,7 @@ def proton_collision_deexcitation_rate(self) -> u.cm**3 / u.s:
494
505
--------
495
506
proton_collision_excitation_rate : Excitation rate due to collisions with protons
496
507
"""
497
- kBTE = np .outer (const . k_B * self .temperature , 1.0 / self ._psplups ['delta_energy' ])
508
+ kBTE = np .outer (self .thermal_energy , 1.0 / self ._psplups ['delta_energy' ])
498
509
omega_upper = 2. * self ._elvlc ['J' ][self ._psplups ['upper_level' ] - 1 ] + 1.
499
510
omega_lower = 2. * self ._elvlc ['J' ][self ._psplups ['lower_level' ] - 1 ] + 1.
500
511
dex_rate = (omega_lower / omega_upper ) * self .proton_collision_excitation_rate * np .exp (1. / kBTE )
@@ -1006,7 +1017,7 @@ def direct_ionization_rate(self) -> u.cm**3 / u.s:
1006
1017
direct_ionization_cross_section : Calculation of :math:`\sigma_I` as a function of :math:`E`.
1007
1018
"""
1008
1019
xgl , wgl = np .polynomial .laguerre .laggauss (12 )
1009
- kBT = const . k_B * self .temperature
1020
+ kBT = self .thermal_energy
1010
1021
energy = np .outer (xgl , kBT ) + self .ionization_potential
1011
1022
cross_section = self .direct_ionization_cross_section (energy )
1012
1023
term1 = np .sqrt (8. / np .pi / const .m_e )* np .sqrt (kBT )* np .exp (- self .ionization_potential / kBT )
@@ -1126,8 +1137,8 @@ def excitation_autoionization_rate(self) -> u.cm**3 / u.s:
1126
1137
Additionally, note that the constant has been rewritten in terms of :math:`h`
1127
1138
rather than :math:`I_H` and :math:`a_0`.
1128
1139
"""
1129
- c = ( const .h ** 2 ) / (( 2. * np .pi * const .m_e )** (1.5 ) * np . sqrt ( const . k_B ) )
1130
- kBTE = np .outer (const . k_B * self .temperature , 1.0 / self ._easplups ['delta_energy' ])
1140
+ c = const .h ** 2 / ( 2. * np .pi * const .m_e )** (1.5 )
1141
+ kBTE = np .outer (self .thermal_energy , 1.0 / self ._easplups ['delta_energy' ])
1131
1142
# NOTE: Transpose here to make final dimensions compatible with multiplication with
1132
1143
# temperature when computing rate
1133
1144
kBTE = kBTE .T
@@ -1139,7 +1150,7 @@ def excitation_autoionization_rate(self) -> u.cm**3 / u.s:
1139
1150
self ._easplups ['bt_type' ])
1140
1151
# NOTE: The 1/omega multiplicity factor is already included in the scaled upsilon
1141
1152
# values provided by CHIANTI
1142
- rate = c * upsilon * np .exp (- 1 / kBTE ) / np .sqrt (self .temperature )
1153
+ rate = c * upsilon * np .exp (- 1 / kBTE ) / np .sqrt (self .thermal_energy )
1143
1154
1144
1155
return rate .sum (axis = 0 )
1145
1156
@@ -1369,13 +1380,13 @@ def free_free(self, wavelength: u.angstrom) -> u.erg * u.cm**3 / u.s / u.angstro
1369
1380
fiasco.IonCollection.free_free: Includes abundance and ionization equilibrium.
1370
1381
"""
1371
1382
prefactor = (const .c / 3. / const .m_e * (const .alpha * const .h / np .pi )** 3
1372
- * np .sqrt (2. * np .pi / 3. / const .m_e / const . k_B ))
1373
- tmp = np .outer (self .temperature , wavelength )
1374
- exp_factor = np .exp (- const .h * const .c / const . k_B / tmp ) / (wavelength ** 2 )
1383
+ * np .sqrt (2. * np .pi / 3. / const .m_e ))
1384
+ tmp = np .outer (self .thermal_energy , wavelength )
1385
+ exp_factor = np .exp (- const .h * const .c / tmp ) / (wavelength ** 2 )
1375
1386
gf = self .gaunt_factor .free_free (self .temperature , wavelength , self .atomic_number , self .charge_state , )
1376
1387
1377
1388
return (prefactor * self .charge_state ** 2 * exp_factor * gf
1378
- / np .sqrt (self .temperature )[:, np .newaxis ])
1389
+ / np .sqrt (self .thermal_energy )[:, np .newaxis ])
1379
1390
1380
1391
@u .quantity_input
1381
1392
def free_free_radiative_loss (self , use_itoh = False ) -> u .erg * u .cm ** 3 / u .s :
@@ -1414,9 +1425,9 @@ def free_free_radiative_loss(self, use_itoh=False) -> u.erg * u.cm**3 / u.s:
1414
1425
--------
1415
1426
fiasco.GauntFactor.free_free_integrated: Calculation of :math:`\langle g_{t,ff}\rangle`.
1416
1427
"""
1417
- prefactor = (16. / 3 ** 1.5 ) * np .sqrt (2. * np .pi * const . k_B / (const .hbar ** 2 * const .m_e ** 3 )) * (const .e .esu ** 6 / const .c ** 3 )
1428
+ prefactor = (16. / 3 ** 1.5 ) * np .sqrt (2. * np .pi / (const .hbar ** 2 * const .m_e ** 3 )) * (const .e .esu ** 6 / const .c ** 3 )
1418
1429
gf = self .gaunt_factor .free_free_integrated (self .temperature , self .charge_state , use_itoh = use_itoh )
1419
- return (prefactor * self .charge_state ** 2 * gf * np .sqrt (self .temperature ))
1430
+ return (prefactor * self .charge_state ** 2 * gf * np .sqrt (self .thermal_energy ))
1420
1431
1421
1432
@needs_dataset ('fblvl' , 'ip' )
1422
1433
@u .quantity_input
@@ -1465,12 +1476,12 @@ def free_bound(self,
1465
1476
:cite:t:`verner_analytic_1995`.
1466
1477
"""
1467
1478
wavelength = np .atleast_1d (wavelength )
1468
- prefactor = ( 2 / np .sqrt (2 * np .pi )/ (const .h * (const .c ** 3 ) * ( const .m_e * const . k_B ) ** (3 / 2 ) ))
1479
+ prefactor = 2 / np .sqrt (2 * np .pi )/ (const .h * (const .c ** 3 ) * const .m_e ** (3 / 2 ))
1469
1480
recombining = self .next_ion ()
1470
1481
omega_0 = recombining ._fblvl ['multiplicity' ][0 ] if recombining ._has_dataset ('fblvl' ) else 1.0
1471
1482
E_photon = const .h * const .c / wavelength
1472
1483
# Precompute this here to avoid repeated outer product calculations
1473
- exp_energy_ratio = np .exp (- np .outer (1 / ( const . k_B * self .temperature ) , E_photon ))
1484
+ exp_energy_ratio = np .exp (- np .outer (1 / self .thermal_energy , E_photon ))
1474
1485
# Fill in observed energies with theoretical energies
1475
1486
E_obs = self ._fblvl ['E_obs' ]* const .h * const .c
1476
1487
E_th = self ._fblvl ['E_th' ]* const .h * const .c
@@ -1499,13 +1510,13 @@ def free_bound(self,
1499
1510
# At these temperatures, the cross-section is 0 anyway so we can just zero
1500
1511
# these terms. Just multiplying by 0 is not sufficient because 0*inf=inf
1501
1512
with np .errstate (over = 'ignore' , invalid = 'ignore' ):
1502
- exp_ip_ratio = np .exp (E_ionize / ( const . k_B * self .temperature ) )
1513
+ exp_ip_ratio = np .exp (E_ionize / self .thermal_energy )
1503
1514
xs_exp_ip_ratio = np .outer (exp_ip_ratio , cross_section )
1504
1515
xs_exp_ip_ratio [:,cross_section == 0.0 * u .cm ** 2 ] = 0.0 * u .cm ** 2
1505
1516
sum_factor += omega * xs_exp_ip_ratio
1506
1517
1507
1518
return (prefactor
1508
- * np .outer (self .temperature ** (- 3 / 2 ), E_photon ** 5 )
1519
+ * np .outer (self .thermal_energy ** (- 3 / 2 ), E_photon ** 5 )
1509
1520
* exp_energy_ratio
1510
1521
* sum_factor / omega_0 )
1511
1522
@@ -1566,8 +1577,8 @@ def free_bound_radiative_loss(self) -> u.erg * u.cm**3 / u.s:
1566
1577
recombined = self .previous_ion ()
1567
1578
if not recombined ._has_dataset ('fblvl' ):
1568
1579
return u .Quantity (np .zeros (self .temperature .shape ) * u .erg * u .cm ** 3 / u .s )
1569
- C_ff = 64 * np .pi / 3.0 * np .sqrt (np .pi / 6. ) * (const .e .esu ** 6 )/ (const .c ** 2 * const .m_e ** 1.5 * const . k_B ** 0.5 )
1570
- prefactor = C_ff * const . k_B * np .sqrt (self .temperature ) / (const .h * const .c )
1580
+ C_ff = 64 * np .pi / 3.0 * np .sqrt (np .pi / 6. ) * (const .e .esu ** 6 )/ (const .c ** 2 * const .m_e ** 1.5 )
1581
+ prefactor = C_ff * np .sqrt (self .thermal_energy ) / (const .h * const .c )
1571
1582
1572
1583
E_obs = recombined ._fblvl ['E_obs' ]* const .h * const .c
1573
1584
E_th = recombined ._fblvl ['E_th' ]* const .h * const .c
@@ -1587,8 +1598,8 @@ def free_bound_radiative_loss(self) -> u.erg * u.cm**3 / u.s:
1587
1598
n0 ,
1588
1599
recombined .ionization_potential ,
1589
1600
ground_state = False )
1590
- term1 = g_fb0 * np .exp (- const .h * const .c / (const . k_B * self .temperature * wvl_n0 ))
1591
- term2 = g_fb1 * np .exp (- const .h * const .c / (const . k_B * self .temperature * wvl_n1 ))
1601
+ term1 = g_fb0 * np .exp (- const .h * const .c / (self .thermal_energy * wvl_n0 ))
1602
+ term2 = g_fb1 * np .exp (- const .h * const .c / (self .thermal_energy * wvl_n1 ))
1592
1603
1593
1604
return prefactor * (term1 + term2 )
1594
1605
0 commit comments