/
collisions.py
1793 lines (1490 loc) · 64 KB
/
collisions.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
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
"""Functions to calculate transport coefficients.
This module includes a number of functions for handling Coulomb collisions
spanning weakly coupled (low density) to strongly coupled (high density)
regimes.
Coulomb collisions
==================
Coulomb collisions are collisions where the interaction force is conveyed
via the electric field, instead of any kind of contact force. They usually
result in relatively small deflections of particle trajectories. However,
given that there are many charged particles in a plasma, one has to take
into account the cumulative effects of many such collisions.
Coulomb logarithms
==================
Please see the documentation for the `Coulomb logarithm <Coulomb_logarithm>`_
for a review of the many ways in which one can define and calculate
that quantity.
Collision rates
===============
The module gathers a few functions helpful for calculating collision
rates between particles. The most general of these is `collision_frequency`,
while if you need average values for a Maxwellian distribution, try
out `collision_rate_electron_ion` and `collision_rate_ion_ion`. These
use `collision_frequency` under the hood.
Macroscopic properties
======================
These include:
* `Spitzer_resistivity`
* `mobility`
* `Knudsen_number`
* `coupling_parameter`
"""
# python modules
from astropy import units as u
import numpy as np
import warnings
# plasmapy modules
from plasmapy import utils
from plasmapy.constants import (c, m_e, k_B, e, eps0, pi, hbar)
from plasmapy import atomic
from plasmapy.physics import parameters
from plasmapy.physics.quantum import (Wigner_Seitz_radius,
thermal_deBroglie_wavelength,
chemical_potential)
from plasmapy.mathematics import Fermi_integral
from plasmapy.utils import check_quantity, _check_relativistic
__all__ = [
"Coulomb_logarithm",
"impact_parameter_perp",
"impact_parameter",
"collision_frequency",
"Coulomb_cross_section",
"fundamental_electron_collision_freq",
"fundamental_ion_collision_freq",
"mean_free_path",
"Spitzer_resistivity",
"mobility",
"Knudsen_number",
"coupling_parameter",
]
@utils.check_quantity(T={"units": u.K, "can_be_negative": False},
n_e={"units": u.m ** -3})
@atomic.particle_input
def Coulomb_logarithm(T,
n_e,
particles: (atomic.Particle, atomic.Particle),
z_mean=np.nan * u.dimensionless_unscaled,
V=np.nan * u.m / u.s,
method="classical"):
r"""
Estimates the Coulomb logarithm.
Parameters
----------
T : ~astropy.units.Quantity
Temperature in units of temperature or energy per particle,
which is assumed to be equal for both the test particle and
the target particle.
n_e : ~astropy.units.Quantity
The electron density in units convertible to per cubic meter.
particles : tuple
A tuple containing string representations of the test particle
(listed first) and the target particle (listed second).
z_mean : ~astropy.units.Quantity, optional
The average ionization (arithmetic mean) for a plasma where the
a macroscopic description is valid. This is used to recover the
average ion density (given the average ionization and electron
density) for calculating the ion sphere radius for non-classical
impact parameters.
V : ~astropy.units.Quantity, optional
The relative velocity between particles. If not provided,
thermal velocity is assumed: :math:`\mu V^2 \sim 2 k_B T`
where `mu` is the reduced mass.
method: str, optional
Selects which theory to use when calculating the Coulomb
logarithm. Defaults to classical method.
Returns
-------
lnLambda : float or numpy.ndarray
An estimate of the Coulomb logarithm that is accurate to
roughly its reciprocal.
Raises
------
ValueError
If the mass or charge of either particle cannot be found, or
any of the inputs contain incorrect values.
UnitConversionError
If the units on any of the inputs are incorrect.
If the n_e, T, or V are not Quantities.
PhysicsError
If the result is smaller than 1.
RelativityError
If the input velocity is same or greater than the speed
of light.
Warns
-----
~astropy.units.UnitsWarning
If units are not provided, SI units are assumed
~plasmapy.utils.RelativityWarning
If the input velocity is greater than 5% of the speed of
light.
Notes
-----
The classical Coulomb logarithm is given by
.. math::
\ln{\Lambda} \equiv \ln\left( \frac{b_{max}}{b_{min}} \right)
where :math:`b_{min}` and :math:`b_{max}` are the inner and outer
impact parameters for Coulomb collisions [1]_.
The outer impact parameter is given by the Debye length:
:math:`b_{min} = \lambda_D` which is a function of electron
temperature and electron density. At distances greater than the
Debye length, electric fields from other particles will be
screened out due to electrons rearranging themselves.
The choice of inner impact parameter is either the distance of closest
approach for a 90 degree Coulomb collision or the thermal deBroglie
wavelength, whichever is larger. This is because Coulomb-style collisions
cannot occur for impact parameters shorter than the deBroglie
wavelength because quantum effects will change the fundamental
nature of the collision [2]_, [3]_.
Errors associated with the classical Coulomb logarithm are of order its
inverse. If the Coulomb logarithm is of order unity, then the
assumptions made in the standard analysis of Coulomb collisions
are invalid.
For dense plasmas where the classical Coulomb logarithm breaks
down there are various extended methods. These can be found
in D.O. Gericke et al's paper, which has a table summarizing
the methods [4]_. The GMS-1 through GMS-6 methods correspond
to the methods found it that table.
It should be noted that GMS-4 thru GMS-6 modify the Coulomb
logarithm to the form:
.. math::
\ln{\Lambda} \equiv 0.5 \ln\left(1 + \frac{b_{max}^2}{b_{min}^2} \right)
This means the Coulomb logarithm will not break down for Lambda < 0,
which occurs for dense, cold plasmas.
Methods
---
Classical
classical Landau-Spitzer approach. Fails for large coupling
parameter where Lambda can become less than zero.
GMS-1
1st method listed in Table 1 of reference [3]
Landau-Spitzer, but with interpolated bmin instead of bmin
selected between deBroglie wavelength and distance of closest
approach. Fails for large coupling
parameter where Lambda can become less than zero.
GMS-2
2nd method listed in Table 1 of reference [3]
Another Landau-Spitzer like approach, but now bmax is also
being interpolated. The interpolation is between the Debye
length and the ion sphere radius, allowing for descriptions
of dilute plasmas. Fails for large coupling
parameter where Lambda can become less than zero.
3rd method listed in Table 1 of reference [3]
classical Landau-Spitzer fails for argument of Coulomb logarithm
Lambda < 0, therefore a clamp is placed at Lambda_min = 2
GMS-4
4th method listed in Table 1 of reference [3]
Spitzer-like extension to Coulomb logarithm by noting that
Coulomb collisions take hyperbolic trajectories. Removes
divergence for small bmin issue in classical Landau-Spitzer
approach, so bmin can be zero. Also doesn't break down as
Lambda < 0 is now impossible, even when coupling parameter is large.
GMS-5
5th method listed in Table 1 of reference [3]
Similar to GMS-4, but setting bmin as distance of closest approach
and bmax interpolated between Debye length and ion sphere radius.
Lambda < 0 impossible.
GMS-6
6th method listed in Table 1 of reference [3]
Similar to GMS-4 and GMS-5, but using interpolation methods
for both bmin and bmax.
Examples
--------
>>> from astropy import units as u
>>> n = 1e19*u.m**-3
>>> T = 1e6*u.K
>>> particles = ('e', 'p')
>>> Coulomb_logarithm(T, n, particles)
14.545527226436974
>>> Coulomb_logarithm(T, n, particles, V=1e6 * u.m / u.s)
11.363478214139432
References
----------
.. [1] Physics of Fully Ionized Gases, L. Spitzer (1962)
.. [2] Francis, F. Chen. Introduction to plasma physics and controlled
fusion 3rd edition. Ch 5 (Springer 2015).
.. [3] Comparison of Coulomb Collision Rates in the Plasma Physics
and Magnetically Confined Fusion Literature, W. Fundamenski and
O.E. Garcia, EFDA–JET–R(07)01
(http://www.euro-fusionscipub.org/wp-content/uploads/2014/11/EFDR07001.pdf)
.. [4] Dense plasma temperature equilibration in the binary collision
approximation. D. O. Gericke et. al. PRE, 65, 036418 (2002).
DOI: 10.1103/PhysRevE.65.036418
"""
# fetching impact min and max impact parameters
bmin, bmax = impact_parameter(T=T,
n_e=n_e,
particles=particles,
z_mean=z_mean,
V=V,
method=method)
if method in ["classical", "GMS-1", "GMS-2"]:
ln_Lambda = np.log(bmax / bmin)
elif method == "GMS-3":
ln_Lambda = np.log(bmax / bmin)
if np.any(ln_Lambda < 2):
if np.isscalar(ln_Lambda.value):
ln_Lambda = 2 * u.dimensionless_unscaled
else:
ln_Lambda[ln_Lambda < 2] = 2 * u.dimensionless_unscaled
elif method in ["GMS-4", "GMS-5", "GMS-6"]:
ln_Lambda = 0.5 * np.log(1 + bmax ** 2 / bmin ** 2)
else:
raise ValueError("Unknown method! Choose from 'classical' and 'GMS-N', N from 1 to 6.")
# applying dimensionless units
ln_Lambda = ln_Lambda.to(u.dimensionless_unscaled).value
if np.any(ln_Lambda < 2) and method in ["classical", "GMS-1", "GMS-2"]:
warnings.warn(f"Coulomb logarithm is {ln_Lambda} and {method} relies on weak coupling.",
utils.CouplingWarning)
elif np.any(ln_Lambda < 4):
warnings.warn(f"Coulomb logarithm is {ln_Lambda}, you might have strong coupling effects",
utils.CouplingWarning)
return ln_Lambda
@atomic.particle_input
def _boilerPlate(T, particles: (atomic.Particle, atomic.Particle), V):
"""
Some boiler plate code for checking if inputs to functions in
collisions.py are good. Also obtains reduced in mass in a
2 particle collision system along with thermal velocity.
"""
# checking temperature is in correct units
T = T.to(u.K, equivalencies=u.temperature_energy())
masses = [p.mass for p in particles]
charges = [np.abs(p.charge) for p in particles]
# obtaining reduced mass of 2 particle collision system
reduced_mass = atomic.reduced_mass(*particles)
# getting thermal velocity of system if no velocity is given
V = _replaceNanVwithThermalV(V, T, reduced_mass)
_check_relativistic(V, 'V')
return T, masses, charges, reduced_mass, V
def _replaceNanVwithThermalV(V, T, m):
"""
Get thermal velocity of system if no velocity is given, for a given mass.
Handles vector checks for V, you must already know that T and m are okay.
"""
if np.any(V == 0):
raise utils.PhysicsError("You cannot have a collision for zero velocity!")
# getting thermal velocity of system if no velocity is given
if V is None:
V = parameters.thermal_speed(T, mass=m)
elif np.any(np.isnan(V)):
if np.isscalar(V.value) and np.isscalar(T.value):
V = parameters.thermal_speed(T, mass=m)
elif np.isscalar(V.value):
V = parameters.thermal_speed(T, mass=m)
elif np.isscalar(T.value):
V = V.copy()
V[np.isnan(V)] = parameters.thermal_speed(T, mass=m)
else:
V = V.copy()
V[np.isnan(V)] = parameters.thermal_speed(T[np.isnan(V)], mass=m)
return V
@check_quantity(T={"units": u.K, "can_be_negative": False})
@atomic.particle_input
def impact_parameter_perp(T,
particles: (atomic.Particle, atomic.Particle),
V=np.nan * u.m / u.s):
r"""Distance of closest approach for a 90 degree Coulomb collision.
Parameters
----------
T : ~astropy.units.Quantity
Temperature in units of temperature or energy per particle,
which is assumed to be equal for both the test particle and
the target particle
particles : tuple
A tuple containing string representations of the test particle
(listed first) and the target particle (listed second)
V : ~astropy.units.Quantity, optional
The relative velocity between particles. If not provided,
thermal velocity is assumed: :math:`\mu V^2 \sim 2 k_B T`
where `mu` is the reduced mass.
Returns
-------
impact_parameter_perp : float or numpy.ndarray
The distance of closest approach for a 90 degree Coulomb collision.
Raises
------
ValueError
If the mass or charge of either particle cannot be found, or
any of the inputs contain incorrect values.
UnitConversionError
If the units on any of the inputs are incorrect
TypeError
If T, or V are not Quantities.
RelativityError
If the input velocity is same or greater than the speed
of light.
Warns
-----
~astropy.units.UnitsWarning
If units are not provided, SI units are assumed
~plasmapy.utils.RelativityWarning
If the input velocity is greater than 5% of the speed of
light.
Notes
-----
The distance of closest approach, impact_parameter_perp, is given by [1]_
.. math::
b_{\perp} = \frac{Z_1 Z_2}{4 \pi \epsilon_0 m v^2}
Examples
--------
>>> from astropy import units as u
>>> T = 1e6*u.K
>>> particles = ('e', 'p')
>>> impact_parameter_perp(T, particles)
<Quantity 8.35505011e-12 m>
References
----------
.. [1] Francis, F. Chen. Introduction to plasma physics and controlled
fusion 3rd edition. Ch 5 (Springer 2015).
"""
# boiler plate checks
T, masses, charges, reduced_mass, V = _boilerPlate(T=T,
particles=particles,
V=V)
# Corresponds to a deflection of 90 degrees, which is valid when
# classical effects dominate.
# !!!Note: an average ionization parameter will have to be
# included here in the future
bPerp = (charges[0] * charges[1] / (4 * pi * eps0 * reduced_mass * V ** 2))
return bPerp.to(u.m)
@check_quantity(T={"units": u.K, "can_be_negative": False},
n_e={"units": u.m ** -3}
)
def impact_parameter(T,
n_e,
particles,
z_mean=np.nan * u.dimensionless_unscaled,
V=np.nan * u.m / u.s,
method="classical"):
r"""Impact parameters for classical and quantum Coulomb collision
Parameters
----------
T : ~astropy.units.Quantity
Temperature in units of temperature or energy per particle,
which is assumed to be equal for both the test particle and
the target particle
n_e : ~astropy.units.Quantity
The electron density in units convertible to per cubic meter.
particles : tuple
A tuple containing string representations of the test particle
(listed first) and the target particle (listed second)
z_mean : ~astropy.units.Quantity, optional
The average ionization (arithmetic mean) for a plasma where the
a macroscopic description is valid. This is used to recover the
average ion density (given the average ionization and electron
density) for calculating the ion sphere radius for non-classical
impact parameters.
V : ~astropy.units.Quantity, optional
The relative velocity between particles. If not provided,
thermal velocity is assumed: :math:`\mu V^2 \sim 2 k_B T`
where `mu` is the reduced mass.
method: str, optional
Selects which theory to use when calculating the Coulomb
logarithm. Defaults to classical method.
Returns
-------
bmin, bmax : tuple of floats
The minimum and maximum impact parameters (distances) for a
Coulomb collision.
Raises
------
ValueError
If the mass or charge of either particle cannot be found, or
any of the inputs contain incorrect values.
UnitConversionError
If the units on any of the inputs are incorrect
TypeError
If the n_e, T, or V are not Quantities.
RelativityError
If the input velocity is same or greater than the speed
of light.
Warns
-----
~astropy.units.UnitsWarning
If units are not provided, SI units are assumed
~plasmapy.utils.RelativityWarning
If the input velocity is greater than 5% of the speed of
light.
Notes
-----
The minimum and maximum impact parameters may be calculated in a
variety of ways. The maximum impact parameter is typically
the Debye length.
For quantum plasmas the maximum impact parameter can be the
quadratic sum of the debye length and ion radius (Wigner_Seitz) [1]_
.. math::
b_{max} = \left(\lambda_{De}^2 + a_i^2\right)^{1/2}
The minimum impact parameter is typically some combination of the
thermal deBroglie wavelength and the distance of closest approach
for a 90 degree Coulomb collision. A quadratic sum is used for
all GMS methods, except for GMS-5, where b_min is simply set to
the distance of closest approach [1]_.
.. math::
b_{min} = \left(\Lambda_{deBroglie}^2 + \rho_{\perp}^2\right)^{1/2}
Examples
--------
>>> from astropy import units as u
>>> n = 1e19*u.m**-3
>>> T = 1e6*u.K
>>> particles = ('e', 'p')
>>> impact_parameter(T, n, particles)
(<Quantity 1.05163088e-11 m>, <Quantity 2.18225522e-05 m>)
>>> impact_parameter(T, n, particles, V=1e6 * u.m / u.s)
(<Quantity 2.53401778e-10 m>, <Quantity 2.18225522e-05 m>)
References
----------
.. [1] Dense plasma temperature equilibration in the binary collision
approximation. D. O. Gericke et. al. PRE, 65, 036418 (2002).
DOI: 10.1103/PhysRevE.65.036418
"""
# boiler plate checks
T, masses, charges, reduced_mass, V = _boilerPlate(T=T,
particles=particles,
V=V)
# catching error where mean charge state is not given for non-classical
# methods that require the ion density
if method == "GMS-2" or method == "GMS-5" or method == "GMS-6":
if np.isnan(z_mean):
raise ValueError("Must provide a z_mean for GMS-2, GMS-5, and "
"GMS-6 methods.")
# Debye length
lambdaDe = parameters.Debye_length(T, n_e)
# deBroglie wavelength
lambdaBroglie = hbar / (2 * reduced_mass * V)
# distance of closest approach in 90 degree Coulomb collision
bPerp = impact_parameter_perp(T=T,
particles=particles,
V=V)
# obtaining minimum and maximum impact parameters depending on which
# method is requested
if method == "classical":
bmax = lambdaDe
# Coulomb-style collisions will not happen for impact parameters
# shorter than either of these two impact parameters, so we choose
# the larger of these two possibilities. That is, between the
# deBroglie wavelength and the distance of closest approach.
# ARRAY NOTES
# T and V should be guaranteed to be same size inputs from _boilerplate
# therefore, lambdaBroglie and bPerp are either both scalar or both array
# if np.isscalar(bPerp.value) and np.isscalar(lambdaBroglie.value): # both scalar
try: # assume both scalar
if bPerp > lambdaBroglie:
bmin = bPerp
else:
bmin = lambdaBroglie
# else: # both lambdaBroglie and bPerp are arrays
except ValueError: # both lambdaBroglie and bPerp are arrays
bmin = lambdaBroglie
bmin[bPerp > lambdaBroglie] = bPerp[bPerp > lambdaBroglie]
elif method == "GMS-1":
# 1st method listed in Table 1 of reference [1]
# This is just another form of the classical Landau-Spitzer
# approach, but bmin is interpolated between the deBroglie
# wavelength and distance of closest approach.
bmax = lambdaDe
bmin = (lambdaBroglie ** 2 + bPerp ** 2) ** (1 / 2)
elif method == "GMS-2":
# 2nd method listed in Table 1 of reference [1]
# Another Landau-Spitzer like approach, but now bmax is also
# being interpolated. The interpolation is between the Debye
# length and the ion sphere radius, allowing for descriptions
# of dilute plasmas.
# Mean ion density.
n_i = n_e / z_mean
# mean ion sphere radius.
ionRadius = Wigner_Seitz_radius(n_i)
bmax = (lambdaDe ** 2 + ionRadius ** 2) ** (1 / 2)
bmin = (lambdaBroglie ** 2 + bPerp ** 2) ** (1 / 2)
elif method == "GMS-3":
# 3rd method listed in Table 1 of reference [1]
# same as GMS-1, but not Lambda has a clamp at Lambda_min = 2
# where Lambda is the argument to the Coulomb logarithm.
bmax = lambdaDe
bmin = (lambdaBroglie ** 2 + bPerp ** 2) ** (1 / 2)
elif method == "GMS-4":
# 4th method listed in Table 1 of reference [1]
bmax = lambdaDe
bmin = (lambdaBroglie ** 2 + bPerp ** 2) ** (1 / 2)
elif method == "GMS-5":
# 5th method listed in Table 1 of reference [1]
# Mean ion density.
n_i = n_e / z_mean
# mean ion sphere radius.
ionRadius = Wigner_Seitz_radius(n_i)
bmax = (lambdaDe ** 2 + ionRadius ** 2) ** (1 / 2)
bmin = bPerp
elif method == "GMS-6":
# 6th method listed in Table 1 of reference [1]
# Mean ion density.
n_i = n_e / z_mean
# mean ion sphere radius.
ionRadius = Wigner_Seitz_radius(n_i)
bmax = (lambdaDe ** 2 + ionRadius ** 2) ** (1 / 2)
bmin = (lambdaBroglie ** 2 + bPerp ** 2) ** (1 / 2)
else:
raise ValueError(f"Method {method} not found!")
# ARRAY NOTES
# it could be that bmin and bmax have different sizes. If Te is a scalar,
# T and V will be scalar from _boilerplate, so bmin will scalar. However
# if n_e is an array, than bmax will be an array. if this is the case,
# do we want to extend the scalar bmin to equal the length of bmax? Sure.
if np.isscalar(bmin.value) and not np.isscalar(bmax.value):
bmin = np.repeat(bmin, len(bmax))
return bmin.to(u.m), bmax.to(u.m)
@check_quantity(T={"units": u.K, "can_be_negative": False},
n={"units": u.m ** -3}
)
def collision_frequency(T,
n,
particles,
z_mean=np.nan * u.dimensionless_unscaled,
V=np.nan * u.m / u.s,
method="classical"):
r"""Collision frequency of particles in a plasma.
Parameters
----------
T : ~astropy.units.Quantity
Temperature in units of temperature.
This should be the electron temperature for electron-electron
and electron-ion collisions, and the ion temperature for
ion-ion collisions.
n : ~astropy.units.Quantity
The density in units convertible to per cubic meter.
This should be the electron density for electron-electron collisions,
and the ion density for electron-ion and ion-ion collisions.
particles : tuple
A tuple containing string representations of the test particle
(listed first) and the target particle (listed second)
z_mean : ~astropy.units.Quantity, optional
The average ionization (arithmetic mean) for a plasma where the
a macroscopic description is valid. This is used to recover the
average ion density (given the average ionization and electron
density) for calculating the ion sphere radius for non-classical
impact parameters.
V : ~astropy.units.Quantity, optional
The relative velocity between particles. If not provided,
thermal velocity is assumed: :math:`\mu V^2 \sim 2 k_B T`
where `mu` is the reduced mass.
method: str, optional
Selects which theory to use when calculating the Coulomb
logarithm. Defaults to classical method.
Returns
-------
freq : float or numpy.ndarray
The collision frequency of particles in a plasma.
Raises
------
ValueError
If the mass or charge of either particle cannot be found, or
any of the inputs contain incorrect values.
UnitConversionError
If the units on any of the inputs are incorrect
TypeError
If the n_e, T, or V are not Quantities.
RelativityError
If the input velocity is same or greater than the speed
of light.
Warns
-----
~astropy.units.UnitsWarning
If units are not provided, SI units are assumed
~plasmapy.utils.RelativityWarning
If the input velocity is greater than 5% of the speed of
light.
Notes
-----
The collision frequency is given by [1]_
.. math::
\nu = n \sigma v \ln{\Lambda}
where n is the particle density, :math:`\sigma` is the collisional
cross-section, :math:`v` is the inter-particle velocity (typically
taken as the thermal velocity), and :math:`\ln{\Lambda}` is the Coulomb
logarithm accounting for small angle collisions.
See eq (2.14) in [2]_.
Examples
--------
>>> from astropy import units as u
>>> n = 1e19*u.m**-3
>>> T = 1e6*u.K
>>> particles = ('e', 'p')
>>> collision_frequency(T, n, particles)
<Quantity 702505.15998601 Hz>
References
----------
.. [1] Francis, F. Chen. Introduction to plasma physics and controlled
fusion 3rd edition. Ch 5 (Springer 2015).
.. [2] http://homepages.cae.wisc.edu/~callen/chap2.pdf
"""
# boiler plate checks
T, masses, charges, reduced_mass, V_r = _boilerPlate(T=T,
particles=particles,
V=V)
# using a more descriptive name for the thermal velocity using
# reduced mass
V_reduced = V_r
if particles[0] in ('e','e-') and particles[1] in ('e','e-'):
# electron-electron collision
# if a velocity was passed, we use that instead of the reduced
# thermal velocity
V = _replaceNanVwithThermalV(V, T, reduced_mass)
# impact parameter for 90 degree collision
bPerp = impact_parameter_perp(T=T,
particles=particles,
V=V_reduced)
print(T, n, particles, z_mean, method)
# Coulomb logarithm
cou_log = Coulomb_logarithm(T,
n,
particles,
z_mean,
V=np.nan * u.m / u.s,
method=method)
elif particles[0] in ('e','e-') or particles[1] in ('e','e-'):
# electron-ion collision
# Need to manually pass electron thermal velocity to obtain
# correct perpendicular collision radius
# we ignore the reduced velocity and use the electron thermal
# velocity instead
V = _replaceNanVwithThermalV(V, T, m_e)
# need to also correct mass in collision radius from reduced
# mass to electron mass
bPerp = impact_parameter_perp(T=T,
particles=particles,
V=V) * reduced_mass / m_e
# Coulomb logarithm
# !!! may also need to correct Coulomb logarithm to be
# electron-electron version !!!
cou_log = Coulomb_logarithm(T,
n,
particles,
z_mean,
V=np.nan * u.m / u.s,
method=method)
else:
# ion-ion collision
# if a velocity was passed, we use that instead of the reduced
# thermal velocity
V = _replaceNanVwithThermalV(V, T, reduced_mass)
bPerp = impact_parameter_perp(T=T,
particles=particles,
V=V)
# Coulomb logarithm
cou_log = Coulomb_logarithm(T,
n,
particles,
z_mean,
V=np.nan * u.m / u.s,
method=method)
# collisional cross section
sigma = Coulomb_cross_section(bPerp)
# collision frequency where Coulomb logarithm accounts for
# small angle collisions, which are more frequent than large
# angle collisions.
freq = n * sigma * V * cou_log
return freq.to(u.Hz)
@check_quantity(impact_param={'units': u.m, 'can_be_negative': False})
def Coulomb_cross_section(impact_param: u.m):
r"""Cross section for a large angle Coulomb collision.
Parameters
----------
impact_param : ~astropy.units.Quantity
Impact parameter for the collision.
Examples
--------
>>> Coulomb_cross_section(7e-10*u.m)
<Quantity 6.1575216e-18 m2>
>>> Coulomb_cross_section(0.5*u.m)
<Quantity 3.14159265 m2>
Notes
-----
The collisional cross-section (see [1]_ for a graphical demonstration)
for a 90 degree Coulomb collision is obtained by
.. math::
\sigma = \pi (2 * \rho_{\perp})^2
where :math:`\rho_{\perp}` is the distance of closest approach for
a 90 degree Coulomb collision. This function is a generalization of that
calculation. Please note that it is not guaranteed to return the correct
results for small angle collisions.
Returns
-------
~astropy.units.Quantity
The Coulomb collision cross section area.
References
----------
.. [1] https://en.wikipedia.org/wiki/Cross_section_(physics)#Collision_among_gas_particles
"""
sigma = np.pi * (2 * impact_param) ** 2
return sigma
@utils.check_quantity(
T_e={'units': u.K, 'can_be_negative': False},
n_e={'units': u.m ** -3, 'can_be_negative': False}
)
def fundamental_electron_collision_freq(T_e,
n_e,
ion_particle,
coulomb_log=None,
V=None,
coulomb_log_method="classical"):
r"""
Average momentum relaxation rate for a slowly flowing Maxwellian distribution of electrons.
[3]_ provides a derivation of this as an average collision frequency between electrons
and ions for a Maxwellian distribution. It is thus a special case of the collision
frequency with an averaging factor, and is on many occasions in transport theory
the most relevant collision frequency that has to be considered. It is heavily
related to diffusion and resistivity in plasmas.
Parameters
----------
T_e : ~astropy.units.Quantity
The electron temperature of the Maxwellian test electrons
n_e : ~astropy.units.Quantity
The number density of the Maxwellian test electrons
ion_particle: str
String signifying a particle type of the field ions, including charge
state information.
V : ~astropy.units.Quantity, optional
The relative velocity between particles. If not provided,
thermal velocity is assumed: :math:`\mu V^2 \sim 2 k_B T`
where `mu` is the reduced mass.
coulomb_log : float or dimensionless ~astropy.units.Quantity, optional
Option to specify a Coulomb logarithm of the electrons on the ions.
If not specified, the Coulomb log will is calculated using the
`~plasmapy.physics.transport.Coulomb_logarithm` function.
coulomb_log_method : string, optional
Method used for Coulomb logarithm calculation (see that function
for more documentation). Choose from "classical" or "GMS-1" to "GMS-6".
Notes
-----
Equations (2.17) and (2.120) in [3]_ provide the original source used
to implement this formula, however, the simplest form that connects our average
collision frequency to the general collision frequency is is this (from 2.17):
.. math::
\nu_e = \frac{4}{3 \sqrt{\pi}} \nu(v_{Te})
Where :math:`\nu` is the general collision frequency and :math:`v_{Te}`
is the electron thermal velocity (the average, for a Maxwellian distribution).
This implementation of the average collision frequency is is equivalent to:
* 1/tau_e from ref [1]_ eqn (2.5e) pp. 215,
* nu_e from ref [2]_ pp. 33,
References
----------
.. [1] Braginskii, S. I. "Transport processes in a plasma." Reviews of
plasma physics 1 (1965): 205.
.. [2] Huba, J. D. "NRL (Naval Research Laboratory) Plasma Formulary,
revised." Naval Research Lab. Report NRL/PU/6790-16-614 (2016).
https://www.nrl.navy.mil/ppd/content/nrl-plasma-formulary
.. [3] J.D. Callen, Fundamentals of Plasma Physics draft material,
Chapter 2, http://homepages.cae.wisc.edu/~callen/chap2.pdf
Examples
--------
>>> from astropy import units as u
>>> fundamental_electron_collision_freq(0.1 * u.eV, 1e6 / u.m ** 3, 'p')
<Quantity 0.00180172 1 / s>
>>> fundamental_electron_collision_freq(1e6 * u.K, 1e6 / u.m ** 3, 'p')
<Quantity 1.07222852e-07 1 / s>
>>> fundamental_electron_collision_freq(100 * u.eV, 1e20 / u.m ** 3, 'p')
<Quantity 3936037.8595928 1 / s>
>>> fundamental_electron_collision_freq(100 * u.eV, 1e20 / u.m ** 3, 'p', coulomb_log_method = 'GMS-1')
<Quantity 3872922.52743562 1 / s>
>>> fundamental_electron_collision_freq(0.1 * u.eV, 1e6 / u.m ** 3, 'p', V = c / 100)
<Quantity 4.41166015e-07 1 / s>
>>> fundamental_electron_collision_freq(100 * u.eV, 1e20 / u.m ** 3, 'p', coulomb_log = 20)
<Quantity 5812633.74935003 1 / s>
See Also
--------
collision_frequency
fundamental_ion_collision_freq
"""
T_e = T_e.to(u.K, equivalencies=u.temperature_energy())
# specify to use electron thermal velocity (most probable), not based on reduced mass
V = _replaceNanVwithThermalV(V, T_e, m_e)
particles = [ion_particle, 'e-']
Z_i = atomic.integer_charge(ion_particle)
nu = collision_frequency(T_e,
n_e,
particles,
z_mean=Z_i,
V=V,
method=coulomb_log_method
)
coeff = 4 / np.sqrt(np.pi) / 3
# accounting for when a Coulomb logarithm value is passed
if np.any(coulomb_log):
cLog = Coulomb_logarithm(T_e,
n_e,
particles,
z_mean=Z_i,
V=np.nan * u.m / u.s,
method=coulomb_log_method)
# dividing out by typical Coulomb logarithm value implicit in
# the collision frequency calculation and replacing with
# the user defined Coulomb logarithm value
nu_mod = nu * coulomb_log / cLog
nu_e = coeff * nu_mod
else:
nu_e = coeff * nu
return nu_e.to(1 / u.s)
@utils.check_quantity(
T_i={'units': u.K, 'can_be_negative': False},
n_i={'units': u.m ** -3, 'can_be_negative': False}
)
def fundamental_ion_collision_freq(T_i,
n_i,
ion_particle,
coulomb_log=None,
V=None,
coulomb_log_method="classical"):
r"""
Average momentum relaxation rate for a slowly flowing Maxwellian distribution of ions.
[3]_ provides a derivation of this as an average collision frequency between ions
and ions for a Maxwellian distribution. It is thus a special case of the collision
frequency with an averaging factor.
Parameters
----------
T_i : ~astropy.units.Quantity