/
streamdf.py
4002 lines (3815 loc) · 154 KB
/
streamdf.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
# The DF of a tidal stream
import copy
import multiprocessing
import warnings
import numpy
import scipy
from packaging.version import parse as parse_version
from scipy import integrate, interpolate, optimize, special
_SCIPY_VERSION = parse_version(scipy.__version__)
if _SCIPY_VERSION < parse_version("0.10"): # pragma: no cover
from scipy.maxentropy import logsumexp
elif _SCIPY_VERSION < parse_version("0.19"): # pragma: no cover
from scipy.misc import logsumexp
else:
from scipy.special import logsumexp
from ..actionAngle.actionAngleIsochroneApprox import dePeriod
from ..orbit import Orbit
from ..potential import flatten as flatten_potential
from ..util import (
ars,
conversion,
coords,
fast_cholesky_invert,
galpyWarning,
multi,
plot,
stable_cho_factor,
)
from ..util._optional_deps import _APY_LOADED, _APY_UNITS
from ..util.conversion import physical_conversion
from .df import df
if _APY_LOADED:
from astropy import units
_INTERPDURINGSETUP = True
_USEINTERP = True
_USESIMPLE = True
# cast a wide net
_TWOPIWRAPS = numpy.arange(-4, 5) * 2.0 * numpy.pi
_labelDict = {
"x": r"$X$",
"y": r"$Y$",
"z": r"$Z$",
"r": r"$R$",
"phi": r"$\phi$",
"vx": r"$V_X$",
"vy": r"$V_Y$",
"vz": r"$V_Z$",
"vr": r"$V_R$",
"vt": r"$V_T$",
"ll": r"$\mathrm{Galactic\ longitude\, (deg)}$",
"bb": r"$\mathrm{Galactic\ latitude\, (deg)}$",
"dist": r"$\mathrm{distance\, (kpc)}$",
"pmll": r"$\mu_l\,(\mathrm{mas\,yr}^{-1})$",
"pmbb": r"$\mu_b\,(\mathrm{mas\,yr}^{-1})$",
"vlos": r"$V_{\mathrm{los}}\,(\mathrm{km\,s}^{-1})$",
}
class streamdf(df):
"""The DF of a tidal stream"""
def __init__(
self,
sigv,
progenitor=None,
pot=None,
aA=None,
useTM=False,
tdisrupt=None,
sigMeanOffset=6.0,
leading=True,
sigangle=None,
deltaAngleTrack=None,
nTrackChunks=None,
nTrackIterations=None,
progIsTrack=False,
ro=None,
vo=None,
Vnorm=None,
Rnorm=None,
R0=8.0,
Zsun=0.0208,
vsun=[-11.1, 8.0 * 30.24, 7.25],
multi=None,
interpTrack=_INTERPDURINGSETUP,
useInterp=_USEINTERP,
nosetup=False,
nospreadsetup=False,
approxConstTrackFreq=False,
useTMHessian=False,
custom_transform=None,
):
"""
Initialize the DF of a tidal stream
Parameters
----------
sigv : float or Quantity
Radial velocity dispersion of the progenitor.
progenitor : galpy.orbit.Orbit
Progenitor orbit as Orbit instance (will be re-integrated, so don't bother integrating the orbit before).
pot : galpy.potential.Potential or list thereof, optional
Potential instance or list thereof.
aA : actionAngle instance
ActionAngle instance used to convert (x,v) to actions. Generally a actionAngleIsochroneApprox instance.
useTM : bool, optional
If set to an actionAngleTorus instance, use this to speed up calculations.
tdisrupt : float or Quantity, optional
Time since start of disruption (default: 5 Gyr).
sigMeanOffset : float, optional
Offset between the mean of the frequencies and the progenitor, in units of the largest eigenvalue of the frequency covariance matrix (along the largest eigenvector), should be positive; to model the trailing part, set leading=False (default: 6.0).
leading : bool, optional
If True, model the leading part of the stream; if False, model the trailing part (default: True).
sigangle : float or Quantity, optional
Estimate of the angle spread of the debris initially (default: sigv/122/[1km/s]=1.8sigv in natural coordinates).
deltaAngleTrack : float or Quantity, optional
Angle to estimate the stream track over (rad; or can be Quantity) (default: None).
nTrackChunks : int, optional
Number of chunks to divide the progenitor track in (default: floor(deltaAngleTrack/0.15)+1).
nTrackIterations : int, optional
Number of iterations to perform when establishing the track; each iteration starts from a previous approximation to the track in (x,v) and calculates a new track based on the deviation between the previous track and the desired track in action-angle coordinates; if not set, an appropriate value is determined based on the magnitude of the misalignment between stream and orbit, with larger numbers of iterations for larger misalignments (default: None).
progIsTrack : bool, optional
If True, then the progenitor (x,v) is actually the (x,v) of the stream track at zero angle separation; useful when initializing with an orbit fit; the progenitor's position will be calculated (default: False).
ro : float or Quantity, optional
Distance scale for translation into internal units (default from configuration file).
vo : float or Quantity, optional
Velocity scale for translation into internal units (default from configuration file).
Vnorm : float or Quantity, optional
Deprecated. Use vo instead (default: None).
Rnorm : float or Quantity, optional
Deprecated. Use ro instead (default: None).
R0 : float or Quantity, optional
Galactocentric radius of the Sun (kpc) (can be different from ro) (default: 8.0).
Zsun : float or Quantity, optional
Sun's height above the plane (kpc) (default: 0.0208).
vsun : numpy.ndarray or Quantity, optional
Sun's motion in cylindrical coordinates (vR positive away from center) (can be Quantity array, but not a list of Quantities) (default: [-11.1, 8.0 * 30.24, 7.25]).
multi : int, optional
If set, use multi-processing (default: None).
interpTrack : bool, optional
Interpolate the stream track while setting up the instance (can be done by hand by calling self._interpolate_stream_track() and self._interpolate_stream_track_aA()) (default: _INTERPDURINGSETUP).
useInterp : bool, optional
Use interpolation by default when calculating approximated frequencies and angles (default: _USEINTERP).
nosetup : bool, optional
If True, don't setup the stream track and anything else that is expensive (default: False).
nospreadsetup : bool, optional
If True, don't setup the spread around the stream track (only for nosetup is False) (default: False).
approxConstTrackFreq : bool, optional
If True, approximate the stream assuming that the frequency is constant along the stream (only works with useTM, for which this leads to a significant speed-up) (default: False).
useTMHessian : bool, optional
If True, compute the basic Hessian dO/dJ_prog using TM; otherwise use aA (default: False).
custom_transform : numpy.ndarray, optional
Matrix implementing the rotation from (ra,dec) to a custom set of sky coordinates (default: None).
Notes
-----
- 2013-09-16 - Started - Bovy (IAS)
- 2013-11-25 - Started over - Bovy (IAS)
"""
if ro is None and not Rnorm is None:
warnings.warn(
"WARNING: Rnorm keyword input to streamdf is deprecated in favor of the standard ro keyword",
galpyWarning,
)
ro = Rnorm
if vo is None and not Vnorm is None:
warnings.warn(
"WARNING: Vnorm keyword input to streamdf is deprecated in favor of the standard vo keyword",
galpyWarning,
)
vo = Vnorm
df.__init__(self, ro=ro, vo=vo)
sigv = conversion.parse_velocity(sigv, vo=self._vo)
self._sigv = sigv
if tdisrupt is None:
self._tdisrupt = 5.0 / conversion.time_in_Gyr(self._vo, self._ro)
else:
self._tdisrupt = conversion.parse_time(tdisrupt, ro=self._ro, vo=self._vo)
self._sigMeanOffset = sigMeanOffset
if pot is None: # pragma: no cover
raise OSError("pot= must be set")
self._pot = flatten_potential(pot)
self._aA = aA
if not self._aA._pot == self._pot:
raise OSError(
"Potential in aA does not appear to be the same as given potential pot"
)
self._check_consistent_units()
if useTM:
self._useTM = True
self._aAT = useTM # confusing, no?
self._approxConstTrackFreq = approxConstTrackFreq
if not self._aAT._pot == self._pot:
raise OSError(
"Potential in useTM=actionAngleTorus instance does not appear to be the same as given potential pot"
)
else:
self._useTM = False
if multi is True: # if set to boolean, enable cpu_count processes
self._multi = multiprocessing.cpu_count()
else:
self._multi = multi
self._progenitor_setup(progenitor, leading, useTMHessian)
sigangle = conversion.parse_angle(sigangle)
deltaAngleTrack = conversion.parse_angle(deltaAngleTrack)
self._offset_setup(sigangle, leading, deltaAngleTrack)
# if progIsTrack, calculate the progenitor that gives a track that is approximately the given orbit
if progIsTrack:
self._setup_progIsTrack()
R0 = conversion.parse_length_kpc(R0)
Zsun = conversion.parse_length_kpc(Zsun)
vsun = conversion.parse_velocity_kms(
numpy.array(vsun) if isinstance(vsun, list) else vsun
)
vsun[0] = conversion.parse_velocity_kms(vsun[0])
vsun[1] = conversion.parse_velocity_kms(vsun[1])
vsun[2] = conversion.parse_velocity_kms(vsun[2])
self._setup_coord_transform(R0, Zsun, vsun, progenitor, custom_transform)
# Determine the stream track
if not nosetup:
self._determine_nTrackIterations(nTrackIterations)
self._determine_stream_track(nTrackChunks)
self._useInterp = useInterp
if interpTrack or self._useInterp:
self._interpolate_stream_track()
self._interpolate_stream_track_aA()
self.calc_stream_lb()
if not nospreadsetup:
self._determine_stream_spread()
return None
def _progenitor_setup(self, progenitor, leading, useTMHessian):
"""The part of the setup relating to the progenitor's orbit"""
# Progenitor orbit: Calculate actions, frequencies, and angles for the progenitor
self._progenitor = progenitor() # call to get new Orbit
# Make sure we do not use physical coordinates
self._progenitor.turn_physical_off()
acfs = self._aA.actionsFreqsAngles(
self._progenitor, _firstFlip=(not leading), use_physical=False
)
self._progenitor_jr = acfs[0][0]
self._progenitor_lz = acfs[1][0]
self._progenitor_jz = acfs[2][0]
self._progenitor_Omegar = acfs[3]
self._progenitor_Omegaphi = acfs[4]
self._progenitor_Omegaz = acfs[5]
self._progenitor_Omega = numpy.array([acfs[3], acfs[4], acfs[5]]).reshape(3)
self._progenitor_angler = acfs[6]
self._progenitor_anglephi = acfs[7]
self._progenitor_anglez = acfs[8]
self._progenitor_angle = numpy.array([acfs[6], acfs[7], acfs[8]]).reshape(3)
# Calculate dO/dJ Jacobian at the progenitor
if useTMHessian:
h, fr, fp, fz, e = self._aAT.hessianFreqs(
self._progenitor_jr, self._progenitor_lz, self._progenitor_jz
)
self._dOdJp = h
# Replace frequencies with TM frequencies
self._progenitor_Omegar = fr
self._progenitor_Omegaphi = fp
self._progenitor_Omegaz = fz
self._progenitor_Omega = numpy.array(
[
self._progenitor_Omegar,
self._progenitor_Omegaphi,
self._progenitor_Omegaz,
]
).reshape(3)
else:
self._dOdJp = calcaAJac(
self._progenitor.vxvv[0], self._aA, dxv=None, dOdJ=True, _initacfs=acfs
)
self._dOdJpInv = numpy.linalg.inv(self._dOdJp)
self._dOdJpEig = numpy.linalg.eig(self._dOdJp)
return None
def _offset_setup(self, sigangle, leading, deltaAngleTrack):
"""The part of the setup related to calculating the stream/progenitor offset"""
# From the progenitor orbit, determine the sigmas in J and angle
self._sigjr = (
(self._progenitor.rap() - self._progenitor.rperi()) / numpy.pi * self._sigv
)
self._siglz = self._progenitor.rperi() * self._sigv
self._sigjz = 2.0 * self._progenitor.zmax() / numpy.pi * self._sigv
# Estimate the frequency covariance matrix from a diagonal J matrix x dOdJ
self._sigjmatrix = numpy.diag(
[self._sigjr**2.0, self._siglz**2.0, self._sigjz**2.0]
)
self._sigomatrix = numpy.dot(
self._dOdJp, numpy.dot(self._sigjmatrix, self._dOdJp.T)
)
# Estimate angle spread as the ratio of the largest to the middle eigenvalue
self._sigomatrixEig = numpy.linalg.eig(self._sigomatrix)
self._sigomatrixEigsortIndx = numpy.argsort(self._sigomatrixEig[0])
self._sortedSigOEig = sorted(self._sigomatrixEig[0])
if sigangle is None:
self._sigangle = self._sigv * 1.8
else:
self._sigangle = sigangle
self._sigangle2 = self._sigangle**2.0
self._lnsigangle = numpy.log(self._sigangle)
# Estimate the frequency mean as lying along the direction of the largest eigenvalue
self._dsigomeanProgDirection = self._sigomatrixEig[1][
:, numpy.argmax(self._sigomatrixEig[0])
]
self._progenitor_Omega_along_dOmega = numpy.dot(
self._progenitor_Omega, self._dsigomeanProgDirection
)
# Make sure we are modeling the correct part of the stream
self._leading = leading
self._sigMeanSign = 1.0
if self._leading and self._progenitor_Omega_along_dOmega < 0.0:
self._sigMeanSign = -1.0
elif not self._leading and self._progenitor_Omega_along_dOmega > 0.0:
self._sigMeanSign = -1.0
self._progenitor_Omega_along_dOmega *= self._sigMeanSign
self._sigomean = (
self._progenitor_Omega
+ self._sigMeanOffset
* self._sigMeanSign
* numpy.sqrt(numpy.amax(self._sigomatrixEig[0]))
* self._dsigomeanProgDirection
)
# numpy.dot(self._dOdJp,
# numpy.array([self._sigjr,self._siglz,self._sigjz]))
self._dsigomeanProg = self._sigomean - self._progenitor_Omega
self._meandO = self._sigMeanOffset * numpy.sqrt(
numpy.amax(self._sigomatrixEig[0])
)
# Store cholesky of sigomatrix for fast evaluation
self._sigomatrixNorm = numpy.sqrt(numpy.sum(self._sigomatrix**2.0))
self._sigomatrixinv, self._sigomatrixLogdet = fast_cholesky_invert(
self._sigomatrix / self._sigomatrixNorm, tiny=10.0**-15.0, logdet=True
)
self._sigomatrixinv /= self._sigomatrixNorm
deltaAngleTrackLim = (
(self._sigMeanOffset + 4.0)
* numpy.sqrt(self._sortedSigOEig[2])
* self._tdisrupt
)
if deltaAngleTrack is None:
deltaAngleTrack = deltaAngleTrackLim
else:
if deltaAngleTrack > deltaAngleTrackLim:
warnings.warn(
"WARNING: angle range large compared to plausible value.",
galpyWarning,
)
self._deltaAngleTrack = deltaAngleTrack
return None
def _setup_coord_transform(self, R0, Zsun, vsun, progenitor, custom_transform):
# Set the coordinate-transformation parameters; check that these do not conflict with those in the progenitor orbit object; need to use the original, since this objects _progenitor has physical turned off
if progenitor._roSet and (
numpy.fabs(self._ro - progenitor._ro) > 10.0**-0.8
or numpy.fabs(R0 - progenitor._ro) > 10.0**-8.0
):
warnings.warn(
"Warning: progenitor's ro does not agree with streamdf's ro and R0; this may have unexpected consequences when projecting into observables",
galpyWarning,
)
if progenitor._voSet and numpy.fabs(self._vo - progenitor._vo) > 10.0**-8.0:
warnings.warn(
"Warning: progenitor's vo does not agree with streamdf's vo; this may have unexpected consequences when projecting into observables",
galpyWarning,
)
if (progenitor._roSet or progenitor._voSet) and numpy.fabs(
Zsun - progenitor._zo
) > 10.0**-8.0:
warnings.warn(
"Warning: progenitor's zo does not agree with streamdf's Zsun; this may have unexpected consequences when projecting into observables",
galpyWarning,
)
if (progenitor._roSet or progenitor._voSet) and numpy.any(
numpy.fabs(
vsun - numpy.array([0.0, self._vo, 0.0]) - progenitor._solarmotion
)
> 10.0**-8.0
):
warnings.warn(
"Warning: progenitor's solarmotion does not agree with streamdf's vsun (after accounting for vo); this may have unexpected consequences when projecting into observables",
galpyWarning,
)
self._R0 = R0
self._Zsun = Zsun
self._vsun = vsun
self._custom_transform = custom_transform
return None
def _setup_progIsTrack(self):
"""If progIsTrack, the progenitor orbit that was passed to the
streamdf initialization is the track at zero angle separation;
this routine computes an actual progenitor position that gives
the desired track given the parameters of the streamdf"""
# We need to flip the sign of the offset, to go to the progenitor
self._sigMeanSign *= -1.0
# Use _determine_stream_track_single to calculate the track-progenitor
# offset at zero angle separation
prog_stream_offset = _determine_stream_track_single(
self._aA,
self._progenitor,
0.0, # time = 0
self._progenitor_angle,
self._sigMeanSign,
self._dsigomeanProgDirection,
lambda x: self.meanOmega(x, use_physical=False),
0.0,
) # angle = 0
# Setup the new progenitor orbit
progenitor = Orbit(prog_stream_offset[3])
# Flip the offset sign again
self._sigMeanSign *= -1.0
# Now re-do the previous setup
self._progenitor_setup(progenitor, self._leading, False)
self._offset_setup(self._sigangle, self._leading, self._deltaAngleTrack)
return None
@physical_conversion("angle", pop=True)
def misalignment(self, isotropic=False, **kwargs):
"""
Calculate the misalignment between the progenitor's frequency
and the direction along which the stream disrupts.
Parameters
----------
isotropic : bool, optional
If True, return the misalignment assuming an isotropic action distribution.
Returns
-------
float
Misalignment in radians.
Warnings
--------
In versions >1.3, the output unit of streamdf.misalignment has been changed to radian (from degree before).
Notes
-----
- 2013-12-05 - Written - Bovy (IAS)
- 2017-10-28 - Changed output unit to rad - Bovy (UofT)
"""
warnings.warn(
"In versions >1.3, the output unit of streamdf.misalignment has been changed to radian (from degree before)",
galpyWarning,
)
if isotropic:
dODir = self._dOdJpEig[1][:, numpy.argmax(numpy.fabs(self._dOdJpEig[0]))]
else:
dODir = self._dsigomeanProgDirection
out = numpy.arccos(
numpy.sum(self._progenitor_Omega * dODir)
/ numpy.sqrt(numpy.sum(self._progenitor_Omega**2.0))
)
if out > numpy.pi / 2.0:
return out - numpy.pi
else:
return out
def freqEigvalRatio(self, isotropic=False):
"""
Calculate the ratio between the largest and 2nd-to-largest (in abs) eigenvalue of sqrt(dO/dJ^T V_J dO/dJ) (if this is big, a 1D stream will form)
Parameters
----------
isotropic : bool, optional
If True, return the ratio assuming an isotropic action distribution (i.e., just of dO/dJ)
Returns
-------
float
Ratio between eigenvalues of fabs(dO / dJ)
Notes
-----
- 2013-12-05 - Written - Bovy (IAS)
"""
if isotropic:
sortedEig = sorted(numpy.fabs(self._dOdJpEig[0]))
return sortedEig[2] / sortedEig[1]
else:
return (
numpy.sqrt(self._sortedSigOEig)[2] / numpy.sqrt(self._sortedSigOEig)[1]
)
@physical_conversion("time", pop=True)
def estimateTdisrupt(self, deltaAngle):
"""
Estimate the time of disruption.
Parameters
----------
deltaAngle : float
Spread in angle since disruption.
Returns
-------
float or Quantity
Time of disruption.
Notes
-----
- 2013-11-27 - Written - Bovy (IAS)
"""
return deltaAngle / numpy.sqrt(numpy.sum(self._dsigomeanProg**2.0))
def subhalo_encounters(
self, venc=numpy.inf, sigma=150.0 / 220.0, nsubhalo=0.3, bmax=0.025, yoon=False
):
"""
Estimate the number of encounters with subhalos over the lifetime of this stream, using a formalism similar to that of Yoon et al. (2011)
Parameters
----------
venc : float, optional
Count encounters with (relative) speeds less than this (relative radial velocity in cylindrical stream frame, unless yoon is True) (can be Quantity)
sigma : float, optional
Velocity dispersion of the DM subhalo population (can be Quantity)
nsubhalo : float, optional
Spatial number density of subhalos (can be Quantity)
bmax : float, optional
Maximum impact parameter (if larger than width of stream) (can be Quantity)
yoon : bool, optional
If True, use erroneous Yoon et al. formula
Returns
-------
float
Number of encounters
Notes
-----
- 2016-01-19 - Written - Bovy (UofT)
"""
venc = conversion.parse_velocity(venc, vo=self._vo)
sigma = conversion.parse_velocity(sigma, vo=self._vo)
nsubhalo = conversion.parse_numdens(nsubhalo, ro=self._ro)
bmax = conversion.parse_length(bmax, ro=self._ro)
Ravg = numpy.mean(
numpy.sqrt(
self._progenitor.orbit[0, :, 0] ** 2.0
+ self._progenitor.orbit[0, :, 3] ** 2.0
)
)
if numpy.isinf(venc):
vencFac = 1.0
elif yoon:
vencFac = 1.0 - (1.0 + venc**2.0 / 4.0 / sigma**2.0) * numpy.exp(
-(venc**2.0) / 4.0 / sigma**2.0
)
else:
vencFac = 1.0 - numpy.exp(-(venc**2.0) / 2.0 / sigma**2.0)
if yoon:
yoonFac = 2 * numpy.sqrt(2.0)
else:
yoonFac = 1.0
# Figure out width of stream
w = self.sigangledAngle(
self._meandO * self._tdisrupt, simple=True, use_physical=False
)
if bmax < w * Ravg / 2.0:
bmax = w * Ravg / 2.0
return (
yoonFac
/ numpy.sqrt(2.0)
* numpy.sqrt(numpy.pi)
* Ravg
* sigma
* self._tdisrupt**2.0
* self._meandO
* bmax
* nsubhalo
* vencFac
)
############################STREAM TRACK FUNCTIONS#############################
def plotTrack(
self, d1="x", d2="z", interp=True, spread=0, simple=_USESIMPLE, *args, **kwargs
):
"""
Plot the stream track.
Parameters
----------
d1 : str, optional
Plot this on the X axis ('x','y','z','R','phi','vx','vy','vz','vR','vt','ll','bb','dist','pmll','pmbb','vlos'). Default is 'x'.
d2 : str, optional
Plot this on the Y axis (same list as for d1). Default is 'z'.
interp : bool, optional
If True, use the interpolated stream track. Default is True.
spread : int, optional
If int > 0, also plot the spread around the track as spread x sigma. Default is 0.
simple : bool, optional
If True, use a simple estimate for the spread in perpendicular angle. Default is _USESIMPLE.
*args : tuple
Arguments passed to galpy.util.plot.plot.
**kwargs : dict
Keyword arguments passed to galpy.util.plot.plot.
Returns
-------
None
Notes
-----
- 2013-12-09 - Written - Bovy (IAS)
"""
if not hasattr(self, "_ObsTrackLB") and (
d1.lower() == "ll"
or d1.lower() == "bb"
or d1.lower() == "dist"
or d1.lower() == "pmll"
or d1.lower() == "pmbb"
or d1.lower() == "vlos"
or d2.lower() == "ll"
or d2.lower() == "bb"
or d2.lower() == "dist"
or d2.lower() == "pmll"
or d2.lower() == "pmbb"
or d2.lower() == "vlos"
):
self.calc_stream_lb()
phys = kwargs.pop("scaleToPhysical", False)
tx = self._parse_track_dim(d1, interp=interp, phys=phys)
ty = self._parse_track_dim(d2, interp=interp, phys=phys)
plot.plot(
tx,
ty,
*args,
xlabel=_labelDict[d1.lower()],
ylabel=_labelDict[d2.lower()],
**kwargs,
)
if spread:
addx, addy = self._parse_track_spread(
d1, d2, interp=interp, phys=phys, simple=simple
)
if ("ls" in kwargs and kwargs["ls"] == "none") or (
"linestyle" in kwargs and kwargs["linestyle"] == "none"
):
kwargs.pop("ls", None)
kwargs.pop("linestyle", None)
spreadls = "none"
else:
spreadls = "-."
spreadmarker = kwargs.pop("marker", None)
spreadcolor = kwargs.pop("color", None)
spreadlw = kwargs.pop("lw", 1.0)
plot.plot(
tx + spread * addx,
ty + spread * addy,
ls=spreadls,
marker=spreadmarker,
color=spreadcolor,
lw=spreadlw,
overplot=True,
)
plot.plot(
tx - spread * addx,
ty - spread * addy,
ls=spreadls,
marker=spreadmarker,
color=spreadcolor,
lw=spreadlw,
overplot=True,
)
return None
def plotProgenitor(self, d1="x", d2="z", *args, **kwargs):
"""
Plot the progenitor orbit.
Parameters
----------
d1 : str, optional
Plot this on the X axis ('x','y','z','R','phi','vx','vy','vz','vR','vt','ll','bb','dist','pmll','pmbb','vlos'). Default is 'x'.
d2 : str, optional
Plot this on the Y axis (same list as for d1). Default is 'z'.
scaleToPhysical : bool, optional
If True, plot positions in kpc and velocities in km/s. Default is False.
*args : tuple
Arguments passed to `galpy.util.plot.plot`.
**kwargs : dict
Keyword arguments passed to `galpy.util.plot.plot`.
Returns
-------
None
Notes
-----
- 2013-12-09 - Written - Bovy (IAS)
"""
tts = self._progenitor.t[
self._progenitor.t < self._trackts[self._nTrackChunks - 1]
]
obs = [self._R0, 0.0, self._Zsun]
obs.extend(self._vsun)
phys = kwargs.pop("scaleToPhysical", False)
tx = self._parse_progenitor_dim(
d1, tts, ro=self._ro, vo=self._vo, obs=obs, phys=phys
)
ty = self._parse_progenitor_dim(
d2, tts, ro=self._ro, vo=self._vo, obs=obs, phys=phys
)
plot.plot(
tx,
ty,
*args,
xlabel=_labelDict[d1.lower()],
ylabel=_labelDict[d2.lower()],
**kwargs,
)
return None
def _parse_track_dim(self, d1, interp=True, phys=False):
"""Parse the dimension to plot the stream track for"""
if interp:
interpStr = "interpolated"
else:
interpStr = ""
if d1.lower() == "x":
tx = self.__dict__["_%sObsTrackXY" % interpStr][:, 0]
elif d1.lower() == "y":
tx = self.__dict__["_%sObsTrackXY" % interpStr][:, 1]
elif d1.lower() == "z":
tx = self.__dict__["_%sObsTrackXY" % interpStr][:, 2]
elif d1.lower() == "r":
tx = self.__dict__["_%sObsTrack" % interpStr][:, 0]
elif d1.lower() == "phi":
tx = self.__dict__["_%sObsTrack" % interpStr][:, 5]
elif d1.lower() == "vx":
tx = self.__dict__["_%sObsTrackXY" % interpStr][:, 3]
elif d1.lower() == "vy":
tx = self.__dict__["_%sObsTrackXY" % interpStr][:, 4]
elif d1.lower() == "vz":
tx = self.__dict__["_%sObsTrackXY" % interpStr][:, 5]
elif d1.lower() == "vr":
tx = self.__dict__["_%sObsTrack" % interpStr][:, 1]
elif d1.lower() == "vt":
tx = self.__dict__["_%sObsTrack" % interpStr][:, 2]
elif d1.lower() == "ll":
tx = self.__dict__["_%sObsTrackLB" % interpStr][:, 0]
elif d1.lower() == "bb":
tx = self.__dict__["_%sObsTrackLB" % interpStr][:, 1]
elif d1.lower() == "dist":
tx = self.__dict__["_%sObsTrackLB" % interpStr][:, 2]
elif d1.lower() == "pmll":
tx = self.__dict__["_%sObsTrackLB" % interpStr][:, 4]
elif d1.lower() == "pmbb":
tx = self.__dict__["_%sObsTrackLB" % interpStr][:, 5]
elif d1.lower() == "vlos":
tx = self.__dict__["_%sObsTrackLB" % interpStr][:, 3]
if phys and (
d1.lower() == "x"
or d1.lower() == "y"
or d1.lower() == "z"
or d1.lower() == "r"
):
tx = copy.copy(tx)
tx *= self._ro
if phys and (
d1.lower() == "vx"
or d1.lower() == "vy"
or d1.lower() == "vz"
or d1.lower() == "vr"
or d1.lower() == "vt"
):
tx = copy.copy(tx)
tx *= self._vo
return tx
def _parse_progenitor_dim(self, d1, ts, ro=None, vo=None, obs=None, phys=False):
"""Parse the dimension to plot the progenitor orbit for"""
if d1.lower() == "x":
tx = self._progenitor.x(ts, ro=ro, vo=vo, obs=obs, use_physical=False)
elif d1.lower() == "y":
tx = self._progenitor.y(ts, ro=ro, vo=vo, obs=obs, use_physical=False)
elif d1.lower() == "z":
tx = self._progenitor.z(ts, ro=ro, vo=vo, obs=obs, use_physical=False)
elif d1.lower() == "r":
tx = self._progenitor.R(ts, ro=ro, vo=vo, obs=obs, use_physical=False)
elif d1.lower() == "phi":
tx = self._progenitor.phi(ts, ro=ro, vo=vo, obs=obs)
elif d1.lower() == "vx":
tx = self._progenitor.vx(ts, ro=ro, vo=vo, obs=obs, use_physical=False)
elif d1.lower() == "vy":
tx = self._progenitor.vy(ts, ro=ro, vo=vo, obs=obs, use_physical=False)
elif d1.lower() == "vz":
tx = self._progenitor.vz(ts, ro=ro, vo=vo, obs=obs, use_physical=False)
elif d1.lower() == "vr":
tx = self._progenitor.vR(ts, ro=ro, vo=vo, obs=obs, use_physical=False)
elif d1.lower() == "vt":
tx = self._progenitor.vT(ts, ro=ro, vo=vo, obs=obs, use_physical=False)
elif d1.lower() == "ll":
tx = self._progenitor.ll(ts, ro=ro, vo=vo, obs=obs)
elif d1.lower() == "bb":
tx = self._progenitor.bb(ts, ro=ro, vo=vo, obs=obs)
elif d1.lower() == "dist":
tx = self._progenitor.dist(ts, ro=ro, vo=vo, obs=obs)
elif d1.lower() == "pmll":
tx = self._progenitor.pmll(ts, ro=ro, vo=vo, obs=obs)
elif d1.lower() == "pmbb":
tx = self._progenitor.pmbb(ts, ro=ro, vo=vo, obs=obs)
elif d1.lower() == "vlos":
tx = self._progenitor.vlos(ts, ro=ro, vo=vo, obs=obs)
if phys and (
d1.lower() == "x"
or d1.lower() == "y"
or d1.lower() == "z"
or d1.lower() == "r"
):
tx = copy.copy(tx)
tx *= self._ro
if phys and (
d1.lower() == "vx"
or d1.lower() == "vy"
or d1.lower() == "vz"
or d1.lower() == "vr"
or d1.lower() == "vt"
):
tx = copy.copy(tx)
tx *= self._vo
return tx
def _parse_track_spread(self, d1, d2, interp=True, phys=False, simple=_USESIMPLE):
"""Determine the spread around the track"""
if not hasattr(self, "_allErrCovs"):
self._determine_stream_spread(simple=simple)
okaySpreadR = ["r", "vr", "vt", "z", "vz", "phi"]
okaySpreadXY = ["x", "y", "z", "vx", "vy", "vz"]
okaySpreadLB = ["ll", "bb", "dist", "vlos", "pmll", "pmbb"]
# Determine which coordinate system we're in
coord = [False, False, False] # R, XY, LB
if d1.lower() in okaySpreadR and d2.lower() in okaySpreadR:
coord[0] = True
elif d1.lower() in okaySpreadXY and d2.lower() in okaySpreadXY:
coord[1] = True
elif d1.lower() in okaySpreadLB and d2.lower() in okaySpreadLB:
coord[2] = True
else:
raise NotImplementedError(
"plotting the spread for coordinates from different systems not implemented yet ..."
)
# Get the right 2D Jacobian
indxDict = {}
indxDict["r"] = 0
indxDict["vr"] = 1
indxDict["vt"] = 2
indxDict["z"] = 3
indxDict["vz"] = 4
indxDict["phi"] = 5
indxDictXY = {}
indxDictXY["x"] = 0
indxDictXY["y"] = 1
indxDictXY["z"] = 2
indxDictXY["vx"] = 3
indxDictXY["vy"] = 4
indxDictXY["vz"] = 5
indxDictLB = {}
indxDictLB["ll"] = 0
indxDictLB["bb"] = 1
indxDictLB["dist"] = 2
indxDictLB["vlos"] = 3
indxDictLB["pmll"] = 4
indxDictLB["pmbb"] = 5
if coord[0]:
relevantCov = self._allErrCovs
relevantDict = indxDict
if phys: # apply scale factors
tcov = copy.copy(relevantCov)
scaleFac = numpy.array(
[self._ro, self._vo, self._vo, self._ro, self._vo, 1.0]
)
tcov *= numpy.tile(scaleFac, (6, 1))
tcov *= numpy.tile(scaleFac, (6, 1)).T
relevantCov = tcov
elif coord[1]:
relevantCov = self._allErrCovsXY
relevantDict = indxDictXY
if phys: # apply scale factors
tcov = copy.copy(relevantCov)
scaleFac = numpy.array(
[self._ro, self._ro, self._ro, self._vo, self._vo, self._vo]
)
tcov *= numpy.tile(scaleFac, (6, 1))
tcov *= numpy.tile(scaleFac, (6, 1)).T
relevantCov = tcov
elif coord[2]:
relevantCov = self._allErrCovsLBUnscaled
relevantDict = indxDictLB
indx0 = numpy.array(
[
[relevantDict[d1.lower()], relevantDict[d1.lower()]],
[relevantDict[d2.lower()], relevantDict[d2.lower()]],
]
)
indx1 = numpy.array(
[
[relevantDict[d1.lower()], relevantDict[d2.lower()]],
[relevantDict[d1.lower()], relevantDict[d2.lower()]],
]
)
cov = relevantCov[:, indx0, indx1] # cov contains all nTrackChunks covs
if not interp:
out = numpy.empty((self._nTrackChunks, 2))
eigDir = numpy.array([1.0, 0.0])
for ii in range(self._nTrackChunks):
covEig = numpy.linalg.eig(cov[ii])
minIndx = numpy.argmin(covEig[0])
minEigvec = covEig[1][
:, minIndx
] # this is the direction of the transverse spread
if numpy.sum(minEigvec * eigDir) < 0.0:
minEigvec *= -1.0 # Keep them pointing in the same direction
out[ii] = minEigvec * numpy.sqrt(covEig[0][minIndx])
eigDir = minEigvec
else:
# We slerp the minor eigenvector and interpolate the eigenvalue
# First store all of the eigenvectors on the track
allEigval = numpy.empty(self._nTrackChunks)
allEigvec = numpy.empty((self._nTrackChunks, 2))
eigDir = numpy.array([1.0, 0.0])
for ii in range(self._nTrackChunks):
covEig = numpy.linalg.eig(cov[ii])
minIndx = numpy.argmin(covEig[0])
minEigvec = covEig[1][
:, minIndx
] # this is the direction of the transverse spread
if numpy.sum(minEigvec * eigDir) < 0.0:
minEigvec *= -1.0 # Keep them pointing in the same direction
allEigval[ii] = numpy.sqrt(covEig[0][minIndx])
allEigvec[ii] = minEigvec
eigDir = minEigvec
# Now interpolate where needed
interpEigval = interpolate.InterpolatedUnivariateSpline(
self._thetasTrack, allEigval, k=3
)
interpolatedEigval = interpEigval(self._interpolatedThetasTrack)
# Interpolate in chunks
interpolatedEigvec = numpy.empty((len(self._interpolatedThetasTrack), 2))
for ii in range(self._nTrackChunks - 1):
slerpOmega = numpy.arccos(numpy.sum(allEigvec[ii] * allEigvec[ii + 1]))
slerpts = (self._interpolatedThetasTrack - self._thetasTrack[ii]) / (
self._thetasTrack[ii + 1] - self._thetasTrack[ii]
)
slerpIndx = (slerpts >= 0.0) * (slerpts <= 1.0)
for jj in range(2):
interpolatedEigvec[slerpIndx, jj] = (
numpy.sin((1 - slerpts[slerpIndx]) * slerpOmega)
* allEigvec[ii, jj]
+ numpy.sin(slerpts[slerpIndx] * slerpOmega)
* allEigvec[ii + 1, jj]
) / numpy.sin(slerpOmega)
out = numpy.tile(interpolatedEigval.T, (2, 1)).T * interpolatedEigvec
if coord[2]: # if LB, undo rescalings that were applied before
out[:, 0] *= self._ErrCovsLBScale[relevantDict[d1.lower()]]
out[:, 1] *= self._ErrCovsLBScale[relevantDict[d2.lower()]]
return (out[:, 0], out[:, 1])
def plotCompareTrackAAModel(self, **kwargs):
"""
Plot the comparison between the underlying model's dOmega_perp vs. dangle_r (line) and the track in (x,v)'s dOmega_perp vs. dangle_r (dots; explicitly calculating the track's action-angle coordinates)
Parameters
----------
**kwargs: dict
galpy.util.plot.plot kwargs
Returns
-------
None
Notes
-----
- 2014-08-27 - Written - Bovy (IAS)
"""
# First calculate the model
model_adiff = (self._ObsTrackAA[:, 3:] - self._progenitor_angle)[
:, 0
] * self._sigMeanSign
model_operp = (
numpy.dot(
self._ObsTrackAA[:, :3] - self._progenitor_Omega,
self._dsigomeanProgDirection,
)
* self._sigMeanSign
)
# Then calculate the track's frequency-angle coordinates
if self._multi is None:
aatrack = numpy.empty((self._nTrackChunks, 6))