-
Notifications
You must be signed in to change notification settings - Fork 306
/
misc.py
298 lines (224 loc) · 7.35 KB
/
misc.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
"""Functions for miscellaneous plasma parameter calculations."""
__all__ = [
"Bohm_diffusion",
"magnetic_energy_density",
"magnetic_pressure",
"thermal_pressure",
]
__aliases__ = ["DB_", "pmag_", "pth_", "ub_"]
import astropy.units as u
from astropy.constants.si import e, k_B, mu0
from plasmapy import particles
from plasmapy.particles import ParticleLike
from plasmapy.utils.decorators import validate_quantities
__all__ += __aliases__
def _grab_charge(ion: ParticleLike, z_mean=None):
"""
Merge two possible inputs for particle charge.
Parameters
----------
ion : `~plasmapy.particles.particle_class.Particle`
a string representing a charged particle, or a Particle object.
z_mean : `float`
An optional float describing the average ionization of a particle
species.
Returns
-------
float
if ``z_mean`` was passed, ``z_mean``, otherwise, the charge number
of ``ion``.
"""
return particles.charge_number(ion) if z_mean is None else z_mean
@validate_quantities(
T_e={"can_be_negative": False, "equivalencies": u.temperature_energy()},
B={"can_be_negative": False},
)
def Bohm_diffusion(
T_e: u.Quantity[u.K], B: u.Quantity[u.T]
) -> u.Quantity[u.m**2 / u.s]:
r"""
Return the Bohm diffusion coefficient.
The Bohm diffusion coefficient was conjectured to follow Bohm model
of the diffusion of plasma across a magnetic field and describe the
diffusion of early fusion energy machines :cite:p:`bohm:1949`. The
rate predicted by Bohm diffusion is much higher than classical
diffusion, and if there were no exceptions, magnetically confined
fusion would be impractical.
.. math::
D_B = \frac{1}{16} \frac{k_B T}{e B}
where :math:`k_B` is the Boltzmann constant
and :math:`e` is the fundamental charge.
**Aliases:** `DB_`
Parameters
----------
T_e : `~astropy.units.Quantity`
The electron temperature.
B : `~astropy.units.Quantity`
The magnitude of the magnetic field in the plasma.
Warns
-----
: `~astropy.units.UnitsWarning`
If units are not provided, SI units are assumed.
Raises
------
`TypeError`
``T_e`` is not a `~astropy.units.Quantity` and cannot be
converted into one.
`~astropy.units.UnitConversionError`
If ``T_e`` is not in appropriate units.
Examples
--------
>>> import astropy.units as u
>>> T_e = 5000 * u.K
>>> B = 10 * u.T
>>> Bohm_diffusion(T_e, B)
<Quantity 0.00269292 m2 / s>
>>> T_e = 50 * u.eV
>>> B = 10 * u.T
>>> Bohm_diffusion(T_e, B)
<Quantity 0.3125 m2 / s>
Returns
-------
D_B : `~astropy.units.Quantity`
The Bohm diffusion coefficient in meters squared per second.
"""
return k_B * T_e / (16 * e * B)
DB_ = Bohm_diffusion
"""Alias to `~plasmapy.formulary.misc.Bohm_diffusion`."""
@validate_quantities
def magnetic_energy_density(B: u.Quantity[u.T]) -> u.Quantity[u.J / u.m**3]:
r"""
Calculate the magnetic energy density.
**Aliases:** `ub_`
Parameters
----------
B : `~astropy.units.Quantity`
The magnetic field in units convertible to tesla.
Returns
-------
E_B : `~astropy.units.Quantity`
The magnetic energy density in units of joules per cubic meter.
Raises
------
`TypeError`
If the input is not a `~astropy.units.Quantity`.
`~astropy.units.UnitConversionError`
If the input is not in units convertible to tesla.
`ValueError`
If the magnetic field strength does not have an appropriate
value.
Warns
-----
: `~astropy.units.UnitsWarning`
If units are not provided, SI units are assumed
Notes
-----
The magnetic energy density is given by:
.. math::
E_B = \frac{B^2}{2 μ_0}
The motivation behind having two separate functions for magnetic
pressure and magnetic energy density is that it allows greater
insight into the physics that are being considered by the user and
thus more readable code.
See Also
--------
magnetic_pressure : Returns an equivalent `~astropy.units.Quantity`,
except in units of pascals.
Examples
--------
>>> import astropy.units as u
>>> magnetic_energy_density(0.1 * u.T)
<Quantity 3978.87... J / m3>
"""
return magnetic_pressure(B)
ub_ = magnetic_energy_density
"""Alias to `~plasmapy.formulary.misc.magnetic_energy_density`."""
@validate_quantities
def magnetic_pressure(B: u.Quantity[u.T]) -> u.Quantity[u.Pa]:
r"""
Calculate the magnetic pressure.
**Aliases:** `pmag_`
Parameters
----------
B : `~astropy.units.Quantity`
The magnetic field in units convertible to tesla.
Returns
-------
p_B : `~astropy.units.Quantity`
The magnetic pressure in units in pascals (newtons per square meter).
Raises
------
`TypeError`
If the input is not a `~astropy.units.Quantity`.
`~astropy.units.UnitConversionError`
If the input is not in units convertible to tesla.
`ValueError`
If the magnetic field strength is not a real number between
:math:`±∞`\ .
Warns
-----
: `~astropy.units.UnitsWarning`
If units are not provided, SI units are assumed.
Notes
-----
The magnetic pressure is given by:
.. math::
p_B = \frac{B^2}{2 μ_0}
The motivation behind having two separate functions for magnetic
pressure and magnetic energy density is that it allows greater
insight into the physics that are being considered by the user and
thus more readable code.
See Also
--------
magnetic_energy_density : returns an equivalent `~astropy.units.Quantity`,
except in units of joules per cubic meter.
Examples
--------
>>> import astropy.units as u
>>> magnetic_pressure(0.1 * u.T).to(u.Pa)
<Quantity 3978.87... Pa>
"""
return (B**2) / (2 * mu0)
pmag_ = magnetic_pressure
"""Alias to `~plasmapy.formulary.misc.magnetic_pressure`."""
@validate_quantities(
T={"can_be_negative": False, "equivalencies": u.temperature_energy()},
n={"can_be_negative": False},
)
def thermal_pressure(T: u.Quantity[u.K], n: u.Quantity[u.m**-3]) -> u.Quantity[u.Pa]:
r"""
Return the thermal pressure for a Maxwellian distribution.
**Aliases:** `pth_`
Parameters
----------
T : `~astropy.units.Quantity`
The particle temperature in either kelvin or energy per particle.
n : `~astropy.units.Quantity`
The particle number density in units convertible to m\ :sup:`-3`\ .
Examples
--------
>>> import astropy.units as u
>>> thermal_pressure(1 * u.eV, 1e20 / u.m**3)
<Quantity 16.021... Pa>
>>> thermal_pressure(10 * u.eV, 1e20 / u.m**3)
<Quantity 160.21... Pa>
Returns
-------
p_th : `~astropy.units.Quantity`
Thermal pressure.
Raises
------
`TypeError`
The temperature or number density is not a `~astropy.units.Quantity`.
`~astropy.units.UnitConversionError`
If the particle temperature is not in units of temperature or
energy per particle.
Notes
-----
The thermal pressure is given by:
.. math::
T_{th} = n k_B T
"""
return n * k_B * T
pth_ = thermal_pressure
"""Alias to `~plasmapy.formulary.misc.thermal_pressure`."""