-
Notifications
You must be signed in to change notification settings - Fork 95
/
actionAngleStaeckel.py
1574 lines (1473 loc) · 53.6 KB
/
actionAngleStaeckel.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
###############################################################################
# actionAngle: a Python module to calculate actions, angles, and frequencies
#
# class: actionAngleStaeckel
#
# Use Binney (2012; MNRAS 426, 1324)'s Staeckel approximation for
# calculating the actions
#
# methods:
# __call__: returns (jr,lz,jz)
#
###############################################################################
import copy
import warnings
import numpy
from scipy import integrate, optimize
from ..potential import (
DiskSCFPotential,
MWPotential,
SCFPotential,
epifreq,
evaluateR2derivs,
evaluateRzderivs,
evaluatez2derivs,
omegac,
verticalfreq,
)
from ..potential.Potential import (
_check_c,
_evaluatePotentials,
_evaluateRforces,
_evaluatezforces,
)
from ..potential.Potential import flatten as flatten_potential
from ..util import coords # for prolate confocal transforms
from ..util import conversion, galpyWarning
from ..util.conversion import physical_conversion, potential_physical_input
from . import actionAngleStaeckel_c
from .actionAngle import UnboundError, actionAngle
from .actionAngleStaeckel_c import _ext_loaded as ext_loaded
class actionAngleStaeckel(actionAngle):
"""Action-angle formalism for axisymmetric potentials using Binney (2012)'s Staeckel approximation"""
def __init__(self, *args, **kwargs):
"""
Initialize an actionAngleStaeckel object.
Parameters
----------
pot : potential or list of potentials (3D)
The potential or list of potentials.
delta : float or Quantity
The focus.
useu0 : bool, optional
Use u0 to calculate dV (not recommended). Default is False.
c : bool, optional
If True, always use C for calculations. Default is False.
order : int, optional
Number of points to use in the Gauss-Legendre numerical integration of the relevant action, frequency, and angle integrals. Default is 10.
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).
Notes
-----
- 2012-11-27 - Started - Bovy (IAS).
"""
actionAngle.__init__(self, ro=kwargs.get("ro", None), vo=kwargs.get("vo", None))
if not "pot" in kwargs: # pragma: no cover
raise OSError("Must specify pot= for actionAngleStaeckel")
self._pot = flatten_potential(kwargs["pot"])
if self._pot == MWPotential:
warnings.warn(
"Use of MWPotential as a Milky-Way-like potential is deprecated; galpy.potential.MWPotential2014, a potential fit to a large variety of dynamical constraints (see Bovy 2015), is the preferred Milky-Way-like potential in galpy",
galpyWarning,
)
if not "delta" in kwargs: # pragma: no cover
raise OSError("Must specify delta= for actionAngleStaeckel")
if ext_loaded and (("c" in kwargs and kwargs["c"]) or not "c" in kwargs):
self._c = _check_c(self._pot)
if "c" in kwargs and kwargs["c"] and not self._c:
warnings.warn(
"C module not used because potential does not have a C implementation",
galpyWarning,
) # pragma: no cover
else:
self._c = False
self._useu0 = kwargs.get("useu0", False)
self._delta = kwargs["delta"]
self._order = kwargs.get("order", 10)
self._delta = conversion.parse_length(self._delta, ro=self._ro)
# Check the units
self._check_consistent_units()
return None
def _evaluate(self, *args, **kwargs):
"""
Evaluate the actions (jr,lz,jz).
Parameters
----------
*args : tuple
Either:
a) R,vR,vT,z,vz[,phi]:
1) floats: phase-space value for single object (phi is optional) (each can be a Quantity)
2) numpy.ndarray: [N] phase-space values for N objects (each can be a Quantity)
b) Orbit instance: initial condition used if that's it, orbit(t) if there is a time given as well as the second argument
delta: bool, optional
can be used to override the object-wide focal length; can also be an array with length N to allow different delta for different phase-space points
u0: float, optional
if object-wide option useu0 is set, u0 to use (if useu0 and useu0 is None, a good value will be computed).
c: bool, optional
True/False to override the object-wide setting for whether or not to use the C implementation.
order: int, optional
number of points to use in the Gauss-Legendre numerical integration of the relevant action integrals.
fixed_quad: bool, optional
if True, use Gaussian quadrature (scipy.integrate.fixed_quad instead of scipy.integrate.quad).
**kwargs: dict, optional
scipy.integrate.fixed_quad or .quad keywords when not using C
Returns
-------
tuple
(jr,lz,jz)
Notes
-----
- 2012-11-27 - Written - Bovy (IAS)
- 2017-12-27 - Allowed individual delta for each point - Bovy (UofT)
"""
delta = kwargs.pop("delta", self._delta)
order = kwargs.get("order", self._order)
if len(args) == 5: # R,vR.vT, z, vz
R, vR, vT, z, vz = args
elif len(args) == 6: # R,vR.vT, z, vz, phi
R, vR, vT, z, vz, phi = args
else:
self._parse_eval_args(*args)
R = self._eval_R
vR = self._eval_vR
vT = self._eval_vT
z = self._eval_z
vz = self._eval_vz
if isinstance(R, float):
R = numpy.array([R])
vR = numpy.array([vR])
vT = numpy.array([vT])
z = numpy.array([z])
vz = numpy.array([vz])
if (
(self._c and not ("c" in kwargs and not kwargs["c"]))
or (ext_loaded and ("c" in kwargs and kwargs["c"]))
) and _check_c(self._pot):
Lz = R * vT
if self._useu0:
# First calculate u0
if "u0" in kwargs:
u0 = numpy.asarray(kwargs["u0"])
else:
E = numpy.array(
[
_evaluatePotentials(self._pot, R[ii], z[ii])
+ vR[ii] ** 2.0 / 2.0
+ vz[ii] ** 2.0 / 2.0
+ vT[ii] ** 2.0 / 2.0
for ii in range(len(R))
]
)
u0 = actionAngleStaeckel_c.actionAngleStaeckel_calcu0(
E, Lz, self._pot, delta
)[0]
kwargs.pop("u0", None)
else:
u0 = None
jr, jz, err = actionAngleStaeckel_c.actionAngleStaeckel_c(
self._pot, delta, R, vR, vT, z, vz, u0=u0, order=order
)
if err == 0:
return (jr, Lz, jz)
else: # pragma: no cover
raise RuntimeError(
"C-code for calculation actions failed; try with c=False"
)
else:
if "c" in kwargs and kwargs["c"] and not self._c: # pragma: no cover
warnings.warn(
"C module not used because potential does not have a C implementation",
galpyWarning,
)
kwargs.pop("c", None)
if len(R) > 1:
ojr = numpy.zeros(len(R))
olz = numpy.zeros(len(R))
ojz = numpy.zeros(len(R))
for ii in range(len(R)):
targs = (R[ii], vR[ii], vT[ii], z[ii], vz[ii])
tkwargs = copy.copy(kwargs)
try:
tkwargs["delta"] = delta[ii]
except (TypeError, IndexError):
tkwargs["delta"] = delta
tjr, tlz, tjz = self(*targs, **tkwargs)
ojr[ii] = tjr[0]
ojz[ii] = tjz[0]
olz[ii] = tlz[0]
return (ojr, olz, ojz)
else:
# Set up the actionAngleStaeckelSingle object
aASingle = actionAngleStaeckelSingle(
R[0], vR[0], vT[0], z[0], vz[0], pot=self._pot, delta=delta
)
return (
numpy.atleast_1d(aASingle.JR(**copy.copy(kwargs))),
numpy.atleast_1d(aASingle._R * aASingle._vT),
numpy.atleast_1d(aASingle.Jz(**copy.copy(kwargs))),
)
def _actionsFreqs(self, *args, **kwargs):
"""
Evaluate the actions and frequencies (jr,lz,jz,Omegar,Omegaphi,Omegaz).
Parameters
----------
*args : tuple
Either:
a) R,vR,vT,z,vz[,phi]:
1) floats: phase-space value for single object (phi is optional) (each can be a Quantity)
2) numpy.ndarray: [N] phase-space values for N objects (each can be a Quantity)
b) Orbit instance: initial condition used if that's it, orbit(t) if there is a time given as well as the second argument
delta: bool, optional
can be used to override the object-wide focal length; can also be an array with length N to allow different delta for different phase-space points
u0: float, optional
if object-wide option useu0 is set, u0 to use (if useu0 and useu0 is None, a good value will be computed).
c: bool, optional
True/False to override the object-wide setting for whether or not to use the C implementation.
order: int, optional
number of points to use in the Gauss-Legendre numerical integration of the relevant action integrals.
fixed_quad: bool, optional
if True, use Gaussian quadrature (scipy.integrate.fixed_quad instead of scipy.integrate.quad).
**kwargs: dict, optional
scipy.integrate.fixed_quad or .quad keywords when not using C
Returns
-------
tuple
(jr,lz,jz,Omegar,Omegaphi,Omegaz)
Notes
-----
- 2013-08-28 - Written - Bovy (IAS)
"""
delta = kwargs.pop("delta", self._delta)
order = kwargs.get("order", self._order)
if (
(self._c and not ("c" in kwargs and not kwargs["c"]))
or (ext_loaded and ("c" in kwargs and kwargs["c"]))
) and _check_c(self._pot):
if len(args) == 5: # R,vR.vT, z, vz
R, vR, vT, z, vz = args
elif len(args) == 6: # R,vR.vT, z, vz, phi
R, vR, vT, z, vz, phi = args
else:
self._parse_eval_args(*args)
R = self._eval_R
vR = self._eval_vR
vT = self._eval_vT
z = self._eval_z
vz = self._eval_vz
if isinstance(R, float):
R = numpy.array([R])
vR = numpy.array([vR])
vT = numpy.array([vT])
z = numpy.array([z])
vz = numpy.array([vz])
Lz = R * vT
if self._useu0:
# First calculate u0
if "u0" in kwargs:
u0 = numpy.asarray(kwargs["u0"])
else:
E = numpy.array(
[
_evaluatePotentials(self._pot, R[ii], z[ii])
+ vR[ii] ** 2.0 / 2.0
+ vz[ii] ** 2.0 / 2.0
+ vT[ii] ** 2.0 / 2.0
for ii in range(len(R))
]
)
u0 = actionAngleStaeckel_c.actionAngleStaeckel_calcu0(
E, Lz, self._pot, delta
)[0]
kwargs.pop("u0", None)
else:
u0 = None
(
jr,
jz,
Omegar,
Omegaphi,
Omegaz,
err,
) = actionAngleStaeckel_c.actionAngleFreqStaeckel_c(
self._pot, delta, R, vR, vT, z, vz, u0=u0, order=order
)
# Adjustments for close-to-circular orbits
indx = numpy.isnan(Omegar) * (jr < 10.0**-3.0) + numpy.isnan(Omegaz) * (
jz < 10.0**-3.0
) # Close-to-circular and close-to-the-plane orbits
if numpy.sum(indx) > 0:
Omegar[indx] = [
epifreq(self._pot, r, use_physical=False) for r in R[indx]
]
Omegaphi[indx] = [
omegac(self._pot, r, use_physical=False) for r in R[indx]
]
Omegaz[indx] = [
verticalfreq(self._pot, r, use_physical=False) for r in R[indx]
]
if err == 0:
return (jr, Lz, jz, Omegar, Omegaphi, Omegaz)
else: # pragma: no cover
raise RuntimeError(
"C-code for calculation actions failed; try with c=False"
)
else:
if "c" in kwargs and kwargs["c"] and not self._c: # pragma: no cover
warnings.warn(
"C module not used because potential does not have a C implementation",
galpyWarning,
)
raise NotImplementedError(
"actionsFreqs with c=False not implemented; maybe you meant to install the C extension?"
)
def _actionsFreqsAngles(self, *args, **kwargs):
"""
Evaluate the actions, frequencies, and angles (jr,lz,jz,Omegar,Omegaphi,Omegaz,angler,anglephi,anglez).
Parameters
----------
*args : tuple
Either:
a) R,vR,vT,z,vz[,phi]:
1) floats: phase-space value for single object (phi is optional) (each can be a Quantity)
2) numpy.ndarray: [N] phase-space values for N objects (each can be a Quantity)
b) Orbit instance: initial condition used if that's it, orbit(t) if there is a time given as well as the second argument
delta: bool, optional
can be used to override the object-wide focal length; can also be an array with length N to allow different delta for different phase-space points
u0: float, optional
if object-wide option useu0 is set, u0 to use (if useu0 and useu0 is None, a good value will be computed).
c: bool, optional
True/False to override the object-wide setting for whether or not to use the C implementation.
order: int, optional
number of points to use in the Gauss-Legendre numerical integration of the relevant action integrals.
fixed_quad: bool, optional
if True, use Gaussian quadrature (scipy.integrate.fixed_quad instead of scipy.integrate.quad).
**kwargs: dict, optional
scipy.integrate.fixed_quad or .quad keywords when not using C
Returns
-------
tuple
(jr,lz,jz,Omegar,Omegaphi,Omegaz,angler,anglephi,anglez)
Notes
-----
- 2013-08-28 - Written - Bovy (IAS)
"""
delta = kwargs.pop("delta", self._delta)
order = kwargs.get("order", self._order)
if (
(self._c and not ("c" in kwargs and not kwargs["c"]))
or (ext_loaded and ("c" in kwargs and kwargs["c"]))
) and _check_c(self._pot):
if len(args) == 5: # R,vR.vT, z, vz pragma: no cover
raise OSError("Must specify phi")
elif len(args) == 6: # R,vR.vT, z, vz, phi
R, vR, vT, z, vz, phi = args
else:
self._parse_eval_args(*args)
R = self._eval_R
vR = self._eval_vR
vT = self._eval_vT
z = self._eval_z
vz = self._eval_vz
phi = self._eval_phi
if isinstance(R, float):
R = numpy.array([R])
vR = numpy.array([vR])
vT = numpy.array([vT])
z = numpy.array([z])
vz = numpy.array([vz])
phi = numpy.array([phi])
Lz = R * vT
if self._useu0:
# First calculate u0
if "u0" in kwargs:
u0 = numpy.asarray(kwargs["u0"])
else:
E = numpy.array(
[
_evaluatePotentials(self._pot, R[ii], z[ii])
+ vR[ii] ** 2.0 / 2.0
+ vz[ii] ** 2.0 / 2.0
+ vT[ii] ** 2.0 / 2.0
for ii in range(len(R))
]
)
u0 = actionAngleStaeckel_c.actionAngleStaeckel_calcu0(
E, Lz, self._pot, delta
)[0]
kwargs.pop("u0", None)
else:
u0 = None
(
jr,
jz,
Omegar,
Omegaphi,
Omegaz,
angler,
anglephi,
anglez,
err,
) = actionAngleStaeckel_c.actionAngleFreqAngleStaeckel_c(
self._pot, delta, R, vR, vT, z, vz, phi, u0=u0, order=order
)
# Adjustments for close-to-circular orbits
indx = numpy.isnan(Omegar) * (jr < 10.0**-3.0) + numpy.isnan(Omegaz) * (
jz < 10.0**-3.0
) # Close-to-circular and close-to-the-plane orbits
if numpy.sum(indx) > 0:
Omegar[indx] = [
epifreq(self._pot, r, use_physical=False) for r in R[indx]
]
Omegaphi[indx] = [
omegac(self._pot, r, use_physical=False) for r in R[indx]
]
Omegaz[indx] = [
verticalfreq(self._pot, r, use_physical=False) for r in R[indx]
]
if err == 0:
return (jr, Lz, jz, Omegar, Omegaphi, Omegaz, angler, anglephi, anglez)
else:
raise RuntimeError(
"C-code for calculation actions failed; try with c=False"
) # pragma: no cover
else: # pragma: no cover
if "c" in kwargs and kwargs["c"] and not self._c: # pragma: no cover
warnings.warn(
"C module not used because potential does not have a C implementation",
galpyWarning,
)
raise NotImplementedError(
"actionsFreqs with c=False not implemented; maybe you meant to install the C extension?"
)
def _EccZmaxRperiRap(self, *args, **kwargs):
"""
Evaluate the eccentricity, maximum height above the plane, peri- and apocenter in the Staeckel approximation.
Parameters
----------
*args : tuple
Either:
a) R,vR,vT,z,vz[,phi]:
1) floats: phase-space value for single object (phi is optional) (each can be a Quantity)
2) numpy.ndarray: [N] phase-space values for N objects (each can be a Quantity)
b) Orbit instance: initial condition used if that's it, orbit(t) if there is a time given as well as the second argument
delta: bool, optional
can be used to override the object-wide focal length; can also be an array with length N to allow different delta for different phase-space points
u0: float, optional
if object-wide option useu0 is set, u0 to use (if useu0 and useu0 is None, a good value will be computed).
c: bool, optional
True/False to override the object-wide setting for whether or not to use the C implementation.
Returns
-------
tuple
(e,zmax,rperi,rap)
Notes
-----
- 2017-12-12 - Written - Bovy (UofT)
"""
delta = kwargs.get("delta", self._delta)
umin, umax, vmin = self._uminumaxvmin(*args, **kwargs)
rperi = coords.uv_to_Rz(umin, numpy.pi / 2.0, delta=delta)[0]
rap_tmp, zmax = coords.uv_to_Rz(umax, vmin, delta=delta)
rap = numpy.sqrt(rap_tmp**2.0 + zmax**2.0)
e = (rap - rperi) / (rap + rperi)
return (e, zmax, rperi, rap)
def _uminumaxvmin(self, *args, **kwargs):
"""
Evaluate u_min, u_max, and v_min in the Staeckel approximation.
Parameters
----------
*args : tuple
Either:
a) R,vR,vT,z,vz[,phi]:
1) floats: phase-space value for single object (phi is optional) (each can be a Quantity)
2) numpy.ndarray: [N] phase-space values for N objects (each can be a Quantity)
b) Orbit instance: initial condition used if that's it, orbit(t) if there is a time given as well as the second argument
delta: bool, optional
can be used to override the object-wide focal length; can also be an array with length N to allow different delta for different phase-space points
u0: float, optional
if object-wide option useu0 is set, u0 to use (if useu0 and useu0 is None, a good value will be computed).
c: bool, optional
True/False to override the object-wide setting for whether or not to use the C implementation.
Returns
-------
tuple
(u_min, u_max, v_min)
Notes
-----
- 2017-12-12 - Written - Bovy (UofT)
"""
delta = numpy.atleast_1d(kwargs.pop("delta", self._delta))
if len(args) == 5: # R,vR.vT, z, vz
R, vR, vT, z, vz = args
elif len(args) == 6: # R,vR.vT, z, vz, phi
R, vR, vT, z, vz, phi = args
else:
self._parse_eval_args(*args)
R = self._eval_R
vR = self._eval_vR
vT = self._eval_vT
z = self._eval_z
vz = self._eval_vz
if isinstance(R, float):
R = numpy.array([R])
vR = numpy.array([vR])
vT = numpy.array([vT])
z = numpy.array([z])
vz = numpy.array([vz])
if (
(self._c and not ("c" in kwargs and not kwargs["c"]))
or (ext_loaded and ("c" in kwargs and kwargs["c"]))
) and _check_c(self._pot):
Lz = R * vT
if self._useu0:
# First calculate u0
if "u0" in kwargs:
u0 = numpy.asarray(kwargs["u0"])
else:
E = numpy.array(
[
_evaluatePotentials(self._pot, R[ii], z[ii])
+ vR[ii] ** 2.0 / 2.0
+ vz[ii] ** 2.0 / 2.0
+ vT[ii] ** 2.0 / 2.0
for ii in range(len(R))
]
)
u0 = actionAngleStaeckel_c.actionAngleStaeckel_calcu0(
E, Lz, self._pot, delta
)[0]
kwargs.pop("u0", None)
else:
u0 = None
(
umin,
umax,
vmin,
err,
) = actionAngleStaeckel_c.actionAngleUminUmaxVminStaeckel_c(
self._pot, delta, R, vR, vT, z, vz, u0=u0
)
if err == 0:
return (umin, umax, vmin)
else: # pragma: no cover
raise RuntimeError(
"C-code for calculation actions failed; try with c=False"
)
else:
if "c" in kwargs and kwargs["c"] and not self._c: # pragma: no cover
warnings.warn(
"C module not used because potential does not have a C implementation",
galpyWarning,
)
kwargs.pop("c", None)
if len(R) > 1:
oumin = numpy.zeros(len(R))
oumax = numpy.zeros(len(R))
ovmin = numpy.zeros(len(R))
for ii in range(len(R)):
targs = (R[ii], vR[ii], vT[ii], z[ii], vz[ii])
tkwargs = copy.copy(kwargs)
tkwargs["delta"] = delta[ii] if len(delta) > 1 else delta[0]
tumin, tumax, tvmin = self._uminumaxvmin(*targs, **tkwargs)
oumin[ii] = tumin[0]
oumax[ii] = tumax[0]
ovmin[ii] = tvmin[0]
return (oumin, oumax, ovmin)
else:
# Set up the actionAngleStaeckelSingle object
aASingle = actionAngleStaeckelSingle(
R[0], vR[0], vT[0], z[0], vz[0], pot=self._pot, delta=delta[0]
)
umin, umax = aASingle.calcUminUmax()
vmin = aASingle.calcVmin()
return (
numpy.atleast_1d(umin),
numpy.atleast_1d(umax),
numpy.atleast_1d(vmin),
)
class actionAngleStaeckelSingle(actionAngle):
"""Action-angle formalism for axisymmetric potentials using Binney (2012)'s Staeckel approximation"""
def __init__(self, *args, **kwargs):
"""
Initialize an actionAngleStaeckelSingle object
Parameters
----------
*args : tuple
Either:
a) R,vR,vT,z,vz[,phi]:
1) floats: phase-space value for single object (phi is optional) (each can be a Quantity)
2) numpy.ndarray: [N] phase-space values for N objects (each can be a Quantity)
b) Orbit instance: initial condition used if that's it, orbit(t) if there is a time given as well as the second argument
pot: Potential or list of Potentials
Potential to use
delta: float, optional
focal length of confocal coordinate system
Notes
-----
- 2012-11-27 - Written - Bovy (IAS)
"""
self._parse_eval_args(*args, _noOrbUnitsCheck=True, **kwargs)
self._R = self._eval_R
self._vR = self._eval_vR
self._vT = self._eval_vT
self._z = self._eval_z
self._vz = self._eval_vz
if not "pot" in kwargs: # pragma: no cover
raise OSError("Must specify pot= for actionAngleStaeckelSingle")
self._pot = kwargs["pot"]
if not "delta" in kwargs: # pragma: no cover
raise OSError("Must specify delta= for actionAngleStaeckel")
self._delta = kwargs["delta"]
# Pre-calculate everything
self._ux, self._vx = coords.Rz_to_uv(self._R, self._z, delta=self._delta)
self._sinvx = numpy.sin(self._vx)
self._cosvx = numpy.cos(self._vx)
self._coshux = numpy.cosh(self._ux)
self._sinhux = numpy.sinh(self._ux)
self._pux = self._delta * (
self._vR * self._coshux * self._sinvx
+ self._vz * self._sinhux * self._cosvx
)
self._pvx = self._delta * (
self._vR * self._sinhux * self._cosvx
- self._vz * self._coshux * self._sinvx
)
EL = self.calcEL()
self._E = EL[0]
self._Lz = EL[1]
# Determine umin and umax
self._u0 = kwargs.pop(
"u0", self._ux
) # u0 as defined by Binney does not matter for a
# single action evaluation, so we don't determine it here
self._sinhu0 = numpy.sinh(self._u0)
self._potu0v0 = potentialStaeckel(self._u0, self._vx, self._pot, self._delta)
self._I3U = (
self._E * self._sinhux**2.0
- self._pux**2.0 / 2.0 / self._delta**2.0
- self._Lz**2.0 / 2.0 / self._delta**2.0 / self._sinhux**2.0
)
self._potupi2 = potentialStaeckel(
self._ux, numpy.pi / 2.0, self._pot, self._delta
)
dV = self._coshux**2.0 * self._potupi2 - (
self._sinhux**2.0 + self._sinvx**2.0
) * potentialStaeckel(self._ux, self._vx, self._pot, self._delta)
self._I3V = (
-self._E * self._sinvx**2.0
+ self._pvx**2.0 / 2.0 / self._delta**2.0
+ self._Lz**2.0 / 2.0 / self._delta**2.0 / self._sinvx**2.0
- dV
)
self.calcUminUmax()
self.calcVmin()
return None
def angleR(self, **kwargs):
raise NotImplementedError(
"'angleR' not yet implemented for Staeckel approximation"
)
def TR(self, **kwargs):
raise NotImplementedError("'TR' not implemented yet for Staeckel approximation")
def Tphi(self, **kwargs):
raise NotImplementedError(
"'Tphi' not implemented yet for Staeckel approxximation"
)
def I(self, **kwargs):
raise NotImplementedError("'I' not implemented yet for Staeckel approxximation")
def Jphi(self): # pragma: no cover
return self._R * self._vT
def JR(self, **kwargs):
"""
Calculate the radial action
Parameters
----------
fixed_quad : bool, optional
If True, use n=10 fixed_quad. Default is False.
**kwargs
scipy.integrate.quad keywords
Returns
-------
float
J_R(R,vT,vT)/ro/vc + estimate of the error (nan for fixed_quad)
Notes
-----
- 2012-11-27 - Written - Bovy (IAS)
"""
if hasattr(self, "_JR"): # pragma: no cover
return self._JR
umin, umax = self.calcUminUmax()
# print self._ux, self._pux, (umax-umin)/umax
if (umax - umin) / umax < 10.0**-6:
return numpy.array([0.0])
order = kwargs.pop("order", 10)
if kwargs.pop("fixed_quad", False):
# factor in next line bc integrand=/2delta^2
self._JR = (
1.0
/ numpy.pi
* numpy.sqrt(2.0)
* self._delta
* integrate.fixed_quad(
_JRStaeckelIntegrand,
umin,
umax,
args=(
self._E,
self._Lz,
self._I3U,
self._delta,
self._u0,
self._sinhu0**2.0,
self._vx,
self._sinvx**2.0,
self._potu0v0,
self._pot,
),
n=order,
**kwargs
)[0]
)
else:
self._JR = (
1.0
/ numpy.pi
* numpy.sqrt(2.0)
* self._delta
* integrate.quad(
_JRStaeckelIntegrand,
umin,
umax,
args=(
self._E,
self._Lz,
self._I3U,
self._delta,
self._u0,
self._sinhu0**2.0,
self._vx,
self._sinvx**2.0,
self._potu0v0,
self._pot,
),
**kwargs
)[0]
)
return self._JR
def Jz(self, **kwargs):
"""
Calculate the vertical action
Parameters
----------
fixed_quad : bool, optional
If True, use n=10 fixed_quad. Default is False.
**kwargs
scipy.integrate.quad keywords
Returns
-------
float
J_z(R,vT,vT)/ro/vc + estimate of the error
Notes
-----
- 2012-11-27 - Written - Bovy (IAS)
"""
if hasattr(self, "_JZ"): # pragma: no cover
return self._JZ
vmin = self.calcVmin()
if (numpy.pi / 2.0 - vmin) < 10.0**-7:
return numpy.array([0.0])
order = kwargs.pop("order", 10)
if kwargs.pop("fixed_quad", False):
# factor in next line bc integrand=/2delta^2
self._JZ = (
2.0
/ numpy.pi
* numpy.sqrt(2.0)
* self._delta
* integrate.fixed_quad(
_JzStaeckelIntegrand,
vmin,
numpy.pi / 2,
args=(
self._E,
self._Lz,
self._I3V,
self._delta,
self._ux,
self._coshux**2.0,
self._sinhux**2.0,
self._potupi2,
self._pot,
),
n=order,
**kwargs
)[0]
)
else:
# factor in next line bc integrand=/2delta^2
self._JZ = (
2.0
/ numpy.pi
* numpy.sqrt(2.0)
* self._delta
* integrate.quad(
_JzStaeckelIntegrand,
vmin,
numpy.pi / 2,
args=(
self._E,
self._Lz,
self._I3V,
self._delta,
self._ux,
self._coshux**2.0,
self._sinhux**2.0,
self._potupi2,
self._pot,
),
**kwargs
)[0]
)
return self._JZ
def calcEL(self, **kwargs):
"""
Calculate the energy and angular momentum.
Parameters
----------
**kwargs : dict
scipy.integrate.quadrature keywords
Returns
-------
tuple
A tuple containing the energy and angular momentum.
Notes
-----
- 2012-11-27 - Written - Bovy (IAS)
"""
E, L = calcELStaeckel(self._R, self._vR, self._vT, self._z, self._vz, self._pot)
return (E, L)
def calcUminUmax(self, **kwargs):
"""
Calculate the u 'apocenter' and 'pericenter'
Returns
-------
tuple
(umin,umax)
Notes
-----
- 2012-11-27 - Written - Bovy (IAS)
"""
if hasattr(self, "_uminumax"): # pragma: no cover
return self._uminumax
E, L = self._E, self._Lz
# Calculate value of the integrand at current point, to check whether
# we are at a turning point
current_val = _JRStaeckelIntegrandSquared(
self._ux,
E,
L,
self._I3U,
self._delta,
self._u0,
self._sinhu0**2.0,
self._vx,
self._sinvx**2.0,
self._potu0v0,
self._pot,
)
if (
numpy.fabs(self._pux) < 1e-7 or numpy.fabs(current_val) < 1e-10
): # We are at umin or umax
eps = 10.0**-8.0
peps = _JRStaeckelIntegrandSquared(
self._ux + eps,
E,
L,
self._I3U,
self._delta,
self._u0,
self._sinhu0**2.0,
self._vx,
self._sinvx**2.0,
self._potu0v0,
self._pot,
)
meps = _JRStaeckelIntegrandSquared(
self._ux - eps,
E,
L,
self._I3U,
self._delta,
self._u0,
self._sinhu0**2.0,
self._vx,
self._sinvx**2.0,
self._potu0v0,
self._pot,
)
if peps < 0.0 and meps > 0.0: # we are at umax
umax = self._ux
rstart, prevr = _uminUmaxFindStart(
self._ux,
E,
L,
self._I3U,
self._delta,
self._u0,
self._sinhu0**2.0,
self._vx,
self._sinvx**2.0,
self._potu0v0,
self._pot,
)
if rstart == 0.0:
umin = 0.0
else:
try:
umin = optimize.brentq(
_JRStaeckelIntegrandSquared,
numpy.atleast_1d(rstart)[0],
numpy.atleast_1d(self._ux)[0] - eps,
(
E,
L,
self._I3U,
self._delta,
self._u0,
self._sinhu0**2.0,
self._vx,
self._sinvx**2.0,
self._potu0v0,
self._pot,
),
maxiter=200,
)
except RuntimeError: # pragma: no cover
raise UnboundError("Orbit seems to be unbound")