1
|
|
|
# coding: utf-8 |
2
|
|
|
# Copyright (c) 2020 Stefan Bender |
3
|
|
|
# |
4
|
|
|
# This file is part of pyeppaurora. |
5
|
|
|
# pyeppaurora is free software: you can redistribute it or modify |
6
|
|
|
# it under the terms of the GNU General Public License as published |
7
|
|
|
# by the Free Software Foundation, version 2. |
8
|
|
|
# See accompanying LICENSE file or http://www.gnu.org/licenses/gpl-2.0.html. |
9
|
|
|
"""Atmospheric ionization rate parametrizations |
10
|
|
|
|
11
|
|
|
Includes the atmospheric ionization rate parametrizations for auroral |
12
|
|
|
and medium-energy electron precipitation, 100 eV--1 MeV [1]_, [2]_, and [3]_. |
13
|
|
|
|
14
|
|
|
.. [1] Roble and Ridley, Ann. Geophys., 5A(6), 369--382, 1987 |
15
|
|
|
.. [2] Fang et al., J. Geophys. Res., 113, A09311, 2008 |
16
|
|
|
.. [3] Fang et al., Geophys. Res. Lett., 37, L22106, 2010 |
17
|
|
|
""" |
18
|
|
|
|
19
|
|
|
import numpy as np |
20
|
|
|
|
21
|
|
|
POLY_F2008 = np.array([ |
22
|
|
|
[ 3.49979e-1, -6.18200e-2, -4.08124e-2, 1.65414e-2], |
23
|
|
|
[ 5.85425e-1, -5.00793e-2, 5.69309e-2, -4.02491e-3], |
24
|
|
|
[ 1.69692e-1, -2.58981e-2, 1.96822e-2, 1.20505e-3], |
25
|
|
|
[-1.22271e-1, -1.15532e-2, 5.37951e-6, 1.20189e-3], |
26
|
|
|
[ 1.57018, 2.87896e-1, -4.14857e-1, 5.18158e-2], |
27
|
|
|
[ 8.83195e-1, 4.31402e-2, -8.33599e-2, 1.02515e-2], |
28
|
|
|
[ 1.90953, -4.74704e-2, -1.80200e-1, 2.46652e-2], |
29
|
|
|
[-1.29566, -2.10952e-1, 2.73106e-1, -2.92752e-2] |
30
|
|
|
]) |
31
|
|
|
|
32
|
|
|
POLY_F2010 = np.array([ |
33
|
|
|
[ 1.24616E+0, 1.45903E+0, -2.42269E-1, 5.95459E-2], |
34
|
|
|
[ 2.23976E+0, -4.22918E-7, 1.36458E-2, 2.53332E-3], |
35
|
|
|
[ 1.41754E+0, 1.44597E-1, 1.70433E-2, 6.39717E-4], |
36
|
|
|
[ 2.48775E-1, -1.50890E-1, 6.30894E-9, 1.23707E-3], |
37
|
|
|
[-4.65119E-1, -1.05081E-1, -8.95701E-2, 1.22450E-2], |
38
|
|
|
[ 3.86019E-1, 1.75430E-3, -7.42960E-4, 4.60881E-4], |
39
|
|
|
[-6.45454E-1, 8.49555E-4, -4.28581E-2, -2.99302E-3], |
40
|
|
|
[ 9.48930E-1, 1.97385E-1, -2.50660E-3, -2.06938E-3] |
41
|
|
|
]) |
42
|
|
|
|
43
|
|
|
vpolyval = np.vectorize(np.polyval, signature='(m,n),()->(n)') |
44
|
|
|
|
45
|
|
|
|
46
|
|
View Code Duplication |
def rr1987(energy, flux, scale_height, rho): |
|
|
|
|
47
|
|
|
"""Atmospheric electron energy dissipation Roble and Ridley, 1987 |
48
|
|
|
|
49
|
|
|
Equations (typo corrected) taken from Fang et al., 2008. |
50
|
|
|
|
51
|
|
|
Parameters |
52
|
|
|
---------- |
53
|
|
|
energy: array_like (M,...) |
54
|
|
|
Characteristic energy E_0 [keV] of the Maxwellian distribution. |
55
|
|
|
flux: array_like (M,...) |
56
|
|
|
Integrated energy flux Q_0 [keV / cm² / s¹] |
57
|
|
|
scale_height: array_like (N,...) |
58
|
|
|
The atmospheric scale heights [cm]. |
59
|
|
|
rho: array_like (N,...) |
60
|
|
|
The atmospheric mass density [g / cm³] |
61
|
|
|
|
62
|
|
|
Returns |
63
|
|
|
------- |
64
|
|
|
en_diss: array_like (M,N) |
65
|
|
|
The dissipated energy profiles. |
66
|
|
|
|
67
|
|
|
References |
68
|
|
|
---------- |
69
|
|
|
|
70
|
|
|
.. [#] Roble and Ridley, Ann. Geophys., 5A(6), 369--382, 1987 |
71
|
|
|
""" |
72
|
|
|
_c1 = 2.11685 |
73
|
|
|
_c2 = 2.97035 |
74
|
|
|
_c3 = 2.09710 |
75
|
|
|
_c4 = 0.74054 |
76
|
|
|
_c5 = 0.58795 |
77
|
|
|
_c6 = 1.72746 |
78
|
|
|
_c7 = 1.37459 |
79
|
|
|
_c8 = 0.93296 |
80
|
|
|
|
81
|
|
|
beta = (rho * scale_height / (4 * 1e-6))**(1 / 1.65) # RR 1987, p. 371 |
82
|
|
|
y = beta / energy # Corrected in Fang et al. 2008 (4) |
83
|
|
|
f_y = (_c1 * (y**_c2) * np.exp(-_c3 * (y**_c4)) + |
84
|
|
|
_c5 * (y**_c6) * np.exp(-_c7 * (y**_c8))) |
85
|
|
|
# Corrected in Fang et al. 2008 (2) |
86
|
|
|
en_diss = 0.5 * flux / scale_height * f_y |
87
|
|
|
return en_diss |
88
|
|
|
|
89
|
|
|
|
90
|
|
View Code Duplication |
def rr1987_mod(energy, flux, scale_height, rho): |
|
|
|
|
91
|
|
|
"""Atmospheric electron energy dissipation Roble and Ridley, 1987 |
92
|
|
|
|
93
|
|
|
Equations (typo corrected) taken from Fang et al., 2008. |
94
|
|
|
Modified polynomial values to get closer to Fang et al., 2008, |
95
|
|
|
origin unknown. |
96
|
|
|
|
97
|
|
|
Parameters |
98
|
|
|
---------- |
99
|
|
|
energy: array_like (M,...) |
100
|
|
|
Characteristic energy E_0 [keV] of the Maxwellian distribution. |
101
|
|
|
flux: array_like (M,...) |
102
|
|
|
Integrated energy flux Q_0 [keV / cm² / s¹] |
103
|
|
|
scale_height: array_like (N,...) |
104
|
|
|
The atmospheric scale heights [cm]. |
105
|
|
|
rho: array_like (N,...) |
106
|
|
|
The atmospheric mass density [g / cm³] |
107
|
|
|
|
108
|
|
|
Returns |
109
|
|
|
------- |
110
|
|
|
en_diss: array_like (M,N) |
111
|
|
|
The dissipated energy profiles. |
112
|
|
|
|
113
|
|
|
References |
114
|
|
|
---------- |
115
|
|
|
|
116
|
|
|
.. [#] Roble and Ridley, Ann. Geophys., 5A(6), 369--382, 1987 |
117
|
|
|
""" |
118
|
|
|
# Modified polynomial, origin unknown |
119
|
|
|
_c1 = 3.233 |
120
|
|
|
_c2 = 2.56588 |
121
|
|
|
_c3 = 2.2541 |
122
|
|
|
_c4 = 0.7297198 |
123
|
|
|
_c5 = 1.106907 |
124
|
|
|
_c6 = 1.71349 |
125
|
|
|
_c7 = 1.8835444 |
126
|
|
|
_c8 = 0.86472135 |
127
|
|
|
|
128
|
|
|
# Fang et al., 2008, Eq. (4) |
129
|
|
|
y = (rho * scale_height / (4.6 * 1e-6))**(1 / 1.65) / energy |
130
|
|
|
f_y = (_c1 * (y**_c2) * np.exp(-_c3 * (y**_c4)) + |
131
|
|
|
_c5 * (y**_c6) * np.exp(-_c7 * (y**_c8))) |
132
|
|
|
# energy dissipated [keV] |
133
|
|
|
en_diss = 0.5 * flux / scale_height * f_y |
134
|
|
|
return en_diss |
135
|
|
|
|
136
|
|
|
|
137
|
|
View Code Duplication |
def fang2008(energy, flux, scale_height, rho, pij=POLY_F2008): |
|
|
|
|
138
|
|
|
"""Atmospheric electron energy dissipation from Fang et al., 2008 |
139
|
|
|
|
140
|
|
|
Ionization profile parametrization as derived in Fang et al., 2008 [#]_. |
141
|
|
|
|
142
|
|
|
Parameters |
143
|
|
|
---------- |
144
|
|
|
energy: array_like (M,...) |
145
|
|
|
Characteristic energy E_0 [keV] of the Maxwellian distribution. |
146
|
|
|
flux: array_like (M,...) |
147
|
|
|
Integrated energy flux Q_0 [keV / cm² / s¹] |
148
|
|
|
scale_height: array_like (N,...) |
149
|
|
|
The atmospheric scale height(s) [cm]. |
150
|
|
|
rho: array_like (N,...) |
151
|
|
|
The atmospheric densities [g / cm³], corresponding to the scale heights. |
152
|
|
|
|
153
|
|
|
Returns |
154
|
|
|
------- |
155
|
|
|
en_diss: array_like (M,N) |
156
|
|
|
The dissipated energy profile(s). |
157
|
|
|
|
158
|
|
|
References |
159
|
|
|
---------- |
160
|
|
|
|
161
|
|
|
.. [#] Fang et al., J. Geophys. Res., 113, A09311, 2008, doi: 10.1029/2008JA013384 |
162
|
|
|
""" |
163
|
|
|
def _f_y(_cc, _y): |
164
|
|
|
# Fang et al., 2008, Eq. (6) |
165
|
|
|
_c = _cc.reshape((8, -1)) |
166
|
|
|
return (_c[0] * (_y**_c[1]) * np.exp(-_c[2] * (_y**_c[3])) + |
167
|
|
|
_c[4] * (_y**_c[5]) * np.exp(-_c[6] * (_y**_c[7]))) |
168
|
|
|
# Fang et al., 2008, Eq. (7) |
169
|
|
|
_cs = np.exp(vpolyval(pij[:, ::-1].T, np.log(energy))).T |
170
|
|
|
# Fang et al., 2008, Eq. (4) |
171
|
|
|
y = (rho * scale_height / (4e-6))**(1 / 1.65) / energy |
172
|
|
|
f_y = _f_y(_cs, y) |
173
|
|
|
# Fang et al., 2008, Eq. (2) |
174
|
|
|
en_diss = 0.5 * f_y * flux / scale_height |
175
|
|
|
return en_diss |
176
|
|
|
|
177
|
|
|
|
178
|
|
View Code Duplication |
def fang2010_mono(energy, flux, scale_height, rho, pij=POLY_F2010): |
|
|
|
|
179
|
|
|
r"""Atmospheric electron energy dissipation from Fang et al., 2010 |
180
|
|
|
|
181
|
|
|
Parametrization for mono-energetic electrons [#]_. |
182
|
|
|
|
183
|
|
|
Parameters |
184
|
|
|
---------- |
185
|
|
|
energy: array_like (M,...) |
186
|
|
|
Characteristic energy E_0 [keV] of the Maxwellian distribution. |
187
|
|
|
flux: array_like (M,...) |
188
|
|
|
Integrated energy flux Q_0 [keV / cm² / s¹] |
189
|
|
|
scale_height: array_like (N,...) |
190
|
|
|
The atmospheric scale heights. |
191
|
|
|
rho: array_like (N,...) |
192
|
|
|
The atmospheric densities, corresponding to the scale heights. |
193
|
|
|
|
194
|
|
|
Returns |
195
|
|
|
------- |
196
|
|
|
en_diss: array_like (M,N) |
197
|
|
|
The dissipated energy profile(s). |
198
|
|
|
|
199
|
|
|
References |
200
|
|
|
---------- |
201
|
|
|
|
202
|
|
|
.. [#] Fang et al., Geophys. Res. Lett., 37, L22106, 2010, doi: 10.1029/2010GL045406 |
203
|
|
|
""" |
204
|
|
|
def _f_y(_cc, _y): |
205
|
|
|
# Fang et al., 2008, Eq. (6), Fang et al., 2010 Eq. (4) |
206
|
|
|
_c = _cc.reshape((8, -1)) |
207
|
|
|
return (_c[0] * (_y**_c[1]) * np.exp(-_c[2] * (_y**_c[3])) + |
208
|
|
|
_c[4] * (_y**_c[5]) * np.exp(-_c[6] * (_y**_c[7]))) |
209
|
|
|
# Fang et al., 2010, Eq. (5) |
210
|
|
|
_cs = np.exp(vpolyval(pij[:, ::-1].T, np.log(energy))).T |
211
|
|
|
# Fang et al., 2010, Eq. (1) |
212
|
|
|
y = 2. / energy * (rho * scale_height / (6e-6))**(0.7) |
213
|
|
|
f_y = _f_y(_cs, y) |
214
|
|
|
# Fang et al., 2008, Eq. (2) |
215
|
|
|
en_diss = f_y * flux / scale_height |
216
|
|
|
return en_diss |
217
|
|
|
|
218
|
|
|
|
219
|
|
View Code Duplication |
def fang2010_spec_int(ens, dfluxes, scale_height, rho, pij=POLY_F2010, axis=-1): |
|
|
|
|
220
|
|
|
r"""Integrate over a given energy spectrum |
221
|
|
|
|
222
|
|
|
Integrates over the mono-energetic parametrization `q` from Fang et al., 2010 |
223
|
|
|
using the given differential particle spectrum `phi`: |
224
|
|
|
|
225
|
|
|
:math:`\int_\text{spec} \phi(E) q(E, Q) E \text{d}E` |
226
|
|
|
|
227
|
|
|
Parameters |
228
|
|
|
---------- |
229
|
|
|
ens: array_like (M,...) |
230
|
|
|
Central (bin) energies of the spectrum |
231
|
|
|
dfluxes: array_like (M,...) |
232
|
|
|
Differential particle fluxes in the given bins |
233
|
|
|
scale_height: array_like (N,...) |
234
|
|
|
The atmospheric scale heights |
235
|
|
|
rho: array_like (N,...) |
236
|
|
|
The atmospheric densities, corresponding to the |
237
|
|
|
scale heights. |
238
|
|
|
|
239
|
|
|
Returns |
240
|
|
|
------- |
241
|
|
|
en_diss: array_like (N) |
242
|
|
|
The dissipated energy profile(s). |
243
|
|
|
""" |
244
|
|
|
ediss_f10 = fang2010_mono( |
245
|
|
|
ens[None, None, :], |
246
|
|
|
dfluxes, |
247
|
|
|
scale_height[..., None], |
248
|
|
|
rho[..., None], |
249
|
|
|
pij=pij, |
250
|
|
|
) |
251
|
|
|
return np.trapz(ediss_f10 * ens, ens, axis=axis) |
252
|
|
|
|
253
|
|
|
|
254
|
|
|
def maxwell_general(en, en_0=10.): |
255
|
|
|
"""Maxwell number flux spectrum as in Fang2008 [1] |
256
|
|
|
|
257
|
|
|
Defined in Fang et al., JGR 2008, Eq. (1). |
258
|
|
|
|
259
|
|
|
Parameters |
260
|
|
|
---------- |
261
|
|
|
en: float |
262
|
|
|
Energy in [keV] |
263
|
|
|
en_0: float, optional |
264
|
|
|
Characteristic energy in [keV], i.e. mode of the distribution. |
265
|
|
|
Default: 10 keV |
266
|
|
|
|
267
|
|
|
Returns |
268
|
|
|
------- |
269
|
|
|
phi: float |
270
|
|
|
Differential hemispherical number flux in [keV-1 cm-2 s-1] |
271
|
|
|
([keV] or scaled by 1 keV-2 cm-2 s-1, e.g. ). |
272
|
|
|
""" |
273
|
|
|
return en * np.exp(-en / en_0) |
274
|
|
|
|
275
|
|
|
|
276
|
|
|
def maxwell_pflux(en, en_0=10.): |
277
|
|
|
"""Maxwell particle flux spectrum as in Fang2008 [1] |
278
|
|
|
|
279
|
|
|
Defined in Fang et al., JGR 2008, Eq. (1). |
280
|
|
|
The total precipitating energy flux is fixed to 1 keV cm-2 s-1, |
281
|
|
|
multiply by Q_0 [keV cm-2 s-1] to scale the particle flux. |
282
|
|
|
|
283
|
|
|
Parameters |
284
|
|
|
---------- |
285
|
|
|
en: float |
286
|
|
|
Energy in [keV] |
287
|
|
|
en_0: float, optional |
288
|
|
|
Characteristic energy in [keV], i.e. mode of the distribution. |
289
|
|
|
Default: 10 keV. |
290
|
|
|
|
291
|
|
|
Returns |
292
|
|
|
------- |
293
|
|
|
phi: float |
294
|
|
|
Hemispherical differential particle flux in [keV-1 cm-2 s-1] |
295
|
|
|
([kev-2] scaled by unit energy flux). |
296
|
|
|
""" |
297
|
|
|
return 0.5 / en_0**3 * maxwell_general(en, en_0) |
298
|
|
|
|
299
|
|
|
|
300
|
|
View Code Duplication |
def fang2010_maxw_int(energy, flux, scale_height, rho, bounds=(0.1, 300.), nstep=128, pij=POLY_F2010): |
|
|
|
|
301
|
|
|
""" |
302
|
|
|
Integrate the mono-energetic parametrization over a Maxwellian |
303
|
|
|
|
304
|
|
|
Parameters |
305
|
|
|
---------- |
306
|
|
|
bounds: tuple, optional |
307
|
|
|
(min, max) [keV] of the integration range to integrate the Maxwellian. |
308
|
|
|
Make sure that this is appropriate to encompass the spectrum. |
309
|
|
|
Default: (0.1, 300.) |
310
|
|
|
nsteps: int, optional |
311
|
|
|
Number of integration steps, default: 128. |
312
|
|
|
""" |
313
|
|
|
bounds_l10 = np.log10(bounds) |
314
|
|
|
ens = np.logspace(*bounds_l10, num=nstep) |
315
|
|
|
dflux = flux * maxwell_pflux(ens[:, None], energy) |
316
|
|
|
return fang2010_spec_int(ens, dflux.T, scale_height, rho, pij=pij, axis=-1) |
317
|
|
|
|