-
Notifications
You must be signed in to change notification settings - Fork 95
/
actionAngleIsochroneInverse.py
222 lines (204 loc) · 8.25 KB
/
actionAngleIsochroneInverse.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
###############################################################################
# actionAngle: a Python module to calculate actions, angles, and frequencies
#
# class: actionAngleIsochroneInverse
#
# Calculate (x,v) coordinates for the Isochrone potential from
# given actions-angle coordinates
#
###############################################################################
import numpy
from scipy import optimize
from ..potential import IsochronePotential
from ..util import conversion
from .actionAngleInverse import actionAngleInverse
class actionAngleIsochroneInverse(actionAngleInverse):
"""Inverse action-angle formalism for the isochrone potential, on the Jphi, Jtheta system of Binney & Tremaine (2008); following McGill & Binney (1990) for transformations"""
def __init__(self, *args, **kwargs):
"""
Initialize an actionAngleIsochroneInverse object.
Parameters
----------
b : float or Quantity, optional
Scale parameter of the isochrone parameter.
ip : galpy.potential.IsochronePotential, optional
Instance of a IsochronePotential.
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
-----
- Either specify b or ip.
- 2017-11-14 - Started - Bovy (UofT)
"""
actionAngleInverse.__init__(self, *args, **kwargs)
if not "b" in kwargs and not "ip" in kwargs: # pragma: no cover
raise OSError("Must specify b= for actionAngleIsochrone")
if "ip" in kwargs:
ip = kwargs["ip"]
if not isinstance(ip, IsochronePotential): # pragma: no cover
raise OSError(
"'Provided ip= does not appear to be an instance of an IsochronePotential"
)
# Check the units
self._pot = ip
self._check_consistent_units()
self.b = ip.b
self.amp = ip._amp
else:
self.b = conversion.parse_length(kwargs["b"], ro=self._ro)
rb = numpy.sqrt(self.b**2.0 + 1.0)
self.amp = (self.b + rb) ** 2.0 * rb
# In case we ever decide to implement this in C...
self._c = False
ext_loaded = False
if ext_loaded and (
("c" in kwargs and kwargs["c"]) or not "c" in kwargs
): # pragma: no cover
self._c = True
else:
self._c = False
if not self._c:
self._ip = IsochronePotential(amp=self.amp, b=self.b)
# Define _pot, because some functions that use actionAngle instances need this
self._pot = IsochronePotential(amp=self.amp, b=self.b)
# Check the units
self._check_consistent_units()
return None
def _evaluate(self, jr, jphi, jz, angler, anglephi, anglez, **kwargs):
"""
Evaluate the phase-space coordinates (x,v) for a number of angles on a single torus.
Parameters
----------
jr : float
Radial action.
jphi : float
Azimuthal action.
jz : float
Vertical action.
angler : numpy.ndarray
Radial angle.
anglephi : numpy.ndarray
Azimuthal angle.
anglez : numpy.ndarray
Vertical angle.
Returns
-------
numpy.ndarray
Phase-space coordinates [R,vR,vT,z,vz,phi].
Notes
-----
- 2017-11-14 - Written - Bovy (UofT).
"""
return self._xvFreqs(jr, jphi, jz, angler, anglephi, anglez, **kwargs)[:6]
def _xvFreqs(self, jr, jphi, jz, angler, anglephi, anglez, **kwargs):
"""
Evaluate the phase-space coordinates (x,v) for a number of angles on a single torus as well as the frequencies.
Parameters
----------
jr : float
Radial action.
jphi : float
Azimuthal action.
jz : float
Vertical action.
angler : numpy.ndarray
Radial angle.
anglephi : numpy.ndarray
Azimuthal angle.
anglez : numpy.ndarray
Vertical angle.
Returns
-------
tuple
A tuple containing the phase-space coordinates (R, vR, vT, z, vz, phi), and the frequencies (OmegaR, Omegaphi, Omegaz).
Notes
-----
- 2017-11-15 - Written - Bovy (UofT).
"""
L = jz + numpy.fabs(jphi) # total angular momentum
L2 = L**2.0
sqrtfourbkL2 = numpy.sqrt(L2 + 4.0 * self.b * self.amp)
H = -2.0 * self.amp**2.0 / (2.0 * jr + L + sqrtfourbkL2) ** 2.0
# Calculate the frequencies
omegar = (-2.0 * H) ** 1.5 / self.amp
omegaz = (1.0 + L / sqrtfourbkL2) / 2.0 * omegar
# Start on getting the coordinates
a = -self.amp / 2.0 / H - self.b
ab = a + self.b
e = numpy.sqrt(1.0 + L2 / (2.0 * H * a**2.0))
# Solve Kepler's-ish equation; ar must be between 0 and 2pi
angler = (numpy.atleast_1d(angler) % (-2.0 * numpy.pi)) % (2.0 * numpy.pi)
anglephi = numpy.atleast_1d(anglephi)
anglez = numpy.atleast_1d(anglez)
eta = numpy.empty(len(angler))
for ii, ar in enumerate(angler):
try:
eta[ii] = optimize.newton(
lambda x: x - a * e / ab * numpy.sin(x) - ar,
0.0,
lambda x: 1 - a * e / ab * numpy.cos(x),
)
except RuntimeError:
# Newton-Raphson did not converge, this has to work,
# bc 0 <= ra < 2pi the following start x have different signs
eta[ii] = optimize.brentq(
lambda x: x - a * e / ab * numpy.sin(x) - ar, 0.0, 2.0 * numpy.pi
)
coseta = numpy.cos(eta)
r = a * numpy.sqrt((1.0 - e * coseta) * (1.0 - e * coseta + 2.0 * self.b / a))
vr = numpy.sqrt(self.amp / ab) * a * e * numpy.sin(eta) / r
taneta2 = numpy.tan(eta / 2.0)
tan11 = numpy.arctan(numpy.sqrt((1.0 + e) / (1.0 - e)) * taneta2)
tan12 = numpy.arctan(
numpy.sqrt((a * (1.0 + e) + 2.0 * self.b) / (a * (1.0 - e) + 2.0 * self.b))
* taneta2
)
tan11[tan11 < 0.0] += numpy.pi
tan12[tan12 < 0.0] += numpy.pi
Lambdaeta = tan11 + L / sqrtfourbkL2 * tan12
psi = anglez - omegaz / omegar * angler + Lambdaeta
lowerl = numpy.sqrt(1.0 - jphi**2.0 / L2)
sintheta = numpy.sin(psi) * lowerl
costheta = numpy.sqrt(1.0 - sintheta**2.0)
vtheta = L * lowerl * numpy.cos(psi) / costheta / r
R = r * costheta
z = r * sintheta
vR = vr * costheta - vtheta * sintheta
vz = vr * sintheta + vtheta * costheta
sinu = sintheta / costheta * jphi / L / lowerl
u = numpy.arcsin(sinu)
u[vtheta < 0.0] = numpy.pi - u[vtheta < 0.0]
phi = anglephi - numpy.sign(jphi) * anglez + u
# For non-inclined orbits, phi == psi
phi[True ^ numpy.isfinite(phi)] = psi[True ^ numpy.isfinite(phi)]
phi = phi % (2.0 * numpy.pi)
phi[phi < 0.0] += 2.0 * numpy.pi
return (R, vR, jphi / R, z, vz, phi, omegar, numpy.sign(jphi) * omegaz, omegaz)
def _Freqs(self, jr, jphi, jz, **kwargs):
"""
Return the frequencies corresponding to a torus
Parameters
----------
jr : float
Radial action
jphi : float
Azimuthal action
jz : float
Vertical action
Returns
-------
tuple
A tuple of three floats representing the frequencies (OmegaR, Omegaphi, Omegaz)
Notes
-----
- 2017-11-15 - Written - Bovy (UofT).
"""
L = jz + numpy.fabs(jphi) # total angular momentum
sqrtfourbkL2 = numpy.sqrt(L**2.0 + 4.0 * self.b * self.amp)
H = -2.0 * self.amp**2.0 / (2.0 * jr + L + sqrtfourbkL2) ** 2.0
# Calculate the frequencies
omegar = (-2.0 * H) ** 1.5 / self.amp
omegaz = (1.0 + L / sqrtfourbkL2) / 2.0 * omegar
return (omegar, numpy.sign(jphi) * omegaz, omegaz)