/
ResonanceEquationsOfMotions.py
1142 lines (958 loc) · 34.1 KB
/
ResonanceEquationsOfMotions.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
import rebound as rb
import reboundx
import numpy as np
import theano
import theano.tensor as T
from exoplanet.theano_ops.kepler import KeplerOp
from warnings import warn
from scipy.optimize import lsq_linear
DEBUG = False
def getOmegaMatrix(n):
"""
Get the 2n x 2n skew-symmetric block matrix:
[0 , I_n]
[-I_n, 0 ]
that appears in Hamilton's equations.
Arguments
---------
n : int
Determines matrix dimension
Returns
-------
numpy.array
"""
return np.vstack(
(
np.concatenate([np.zeros((n,n)),np.eye(n)]).T,
np.concatenate([-np.eye(n),np.zeros((n,n))]).T
)
)
def calc_DisturbingFunction_with_sinf_cosf(alpha,e1,e2,w1,w2,sinf1,cosf1,sinf2,cosf2):
"""
Compute the value of the disturbing function
.. math::
\frac{a'}{|r-r'|} - a'\frac{r.r'}{|r'^3|}
from a set of input orbital elements for coplanar planets.
Arguments
---------
alpha : float
semi-major axis ratio
e1 : float
inner eccentricity
e2 : float
outer eccentricity
w1 : float
inner long. of peri
w2 : float
outer long. of peri
sinf1 : float
sine of inner planet true anomaly
cosf1 : float
cosine of inner planet true anomaly
sinf2 : float
sine of outer planet true anomaly
cosf2 : float
cosine of outer planet true anomaly
Returns
-------
float :
Disturbing function value
"""
r1 = alpha * (1-e1*e1) /(1 + e1 * cosf1)
_x1 = r1 * cosf1
_y1 = r1 * sinf1
Cw1 = T.cos(w1)
Sw1 = T.sin(w1)
x1 = Cw1 * _x1 - Sw1 * _y1
y1 = Sw1 * _x1 + Cw1 * _y1
r2 = (1-e2*e2) /(1 + e2 * cosf2)
_x2 = r2 * cosf2
_y2 = r2 * sinf2
Cw2 = T.cos(w2)
Sw2 = T.sin(w2)
x2 = Cw2 * _x2 - Sw2 * _y2
y2 = Sw2 * _x2 + Cw2 * _y2
# direct term
dx = (x2 - x1)
dy = (y2 - y1)
dr2 = dx*dx + dy*dy
direct = 1 / T.sqrt(dr2)
# indirect term
r1dotr = (x2 * x1 + y2 * y1)
r1sq = x2*x2 + y2*y2
r1_3 = 1 / r1sq / T.sqrt(r1sq)
indirect = -1 * r1dotr * r1_3
return direct+indirect
def get_compiled_Hkep_Hpert_full():
# resonance j and k
j,k = T.lscalars('jk')
s = (j-k) / k
# Planet masses: m1,m2
m1,m2 = T.dscalars(2)
Mstar = 1
eta0 = Mstar
eta1 = eta0+m1
eta2 =eta1+m2
mtilde1 = m1 * (eta0/eta1)
Mtilde1 = Mstar * (eta1/eta0)
mtilde2 = m2 * (eta1/eta2)
Mtilde2 = Mstar * (eta2/eta1)
eps = m1 * m2 / (mtilde1 + mtilde2) / Mstar
beta1 = mtilde1 / (mtilde1 + mtilde2)
beta2 = mtilde2 / (mtilde1 + mtilde2)
gamma = mtilde2/mtilde1
# Dynamical variables:
dyvars = T.vector()
Q,sigma1, sigma2, I1, I2, amd = [dyvars[i] for i in range(6)]
# Set lambda2=0
l2 = T.constant(0.)
l1 = -1 * k * Q
w1 = (1+s) * l2 - s * l1 - sigma1
w2 = (1+s) * l2 - s * l1 - sigma2
Gamma1 = I1
Gamma2 = I2
# Resonant semi-major axis ratio
alpha_res = ((j-k)/j)**(2/3) * ((Mstar + m1) / (Mstar+m2))**(1/3)
P0 = k * ( beta2 - beta1 * T.sqrt(alpha_res) ) / 2
P = P0 - k * (s+1/2) * amd
Ltot = beta1 * T.sqrt(alpha_res) + beta2 - amd
L1 = Ltot/2 - P / k - s * (I1 + I2)
L2 = Ltot/2 + P / k + (1 + s) * (I1 + I2)
a1 = (L1 / beta1 )**2 * eta0 / eta1
e1 = T.sqrt(1-(1-(Gamma1 / L1))**2)
a2 = (L2 / beta2 )**2 * eta1 / eta2
e2 = T.sqrt(1-(1-(Gamma2 / L2))**2)
Hkep = - eta1 * beta1 / (2 * a1) / eta0 - eta2 * beta2 / (2 * a2) / eta1
alpha = a1 / a2
ko = KeplerOp()
M1 = l1 - w1
M2 = l2 - w2
sinf1,cosf1 = ko( M1, e1 + T.zeros_like(M1) )
sinf2,cosf2 = ko( M2, e2 + T.zeros_like(M2) )
R = calc_DisturbingFunction_with_sinf_cosf(alpha,e1,e2,w1,w2,sinf1,cosf1,sinf2,cosf2)
Hpert = -eps * R / a2
omega_syn = T.sqrt(eta2/eta1)/a2**1.5 - T.sqrt(eta1/eta0) / a1**1.5
gradHpert = T.grad(Hpert,wrt=dyvars)
gradHkep = T.grad(Hkep,wrt=dyvars)
grad_omega_syn = T.grad(omega_syn,wrt=dyvars)
extra_ins = [m1,m2,j,k]
ins = [dyvars] + extra_ins
# Scalars
omega_syn_fn = theano.function(
inputs=ins,
outputs=omega_syn,
givens=None,
on_unused_input='ignore'
)
Hpert_fn = theano.function(
inputs=ins,
outputs=Hpert,
givens=None,
on_unused_input='ignore'
)
Hkep_fn = theano.function(
inputs=ins,
outputs=Hkep,
givens=None,
on_unused_input='ignore'
)
# gradients
grad_omega_syn_fn = theano.function(
inputs=ins,
outputs=grad_omega_syn,
givens=None,
on_unused_input='ignore'
)
gradHpert_fn = theano.function(
inputs=ins,
outputs=gradHpert,
givens=None,
on_unused_input='ignore'
)
gradHkep_fn = theano.function(
inputs=ins,
outputs=gradHkep,
givens=None,
on_unused_input='ignore'
)
return omega_syn_fn,Hkep_fn,Hpert_fn,grad_omega_syn_fn,gradHkep_fn,gradHpert_fn
def get_compiled_theano_functions(N_QUAD_PTS):
# resonance j and k
j,k = T.lscalars('jk')
s = (j-k) / k
# Planet masses: m1,m2
m1,m2 = T.dscalars(2)
Mstar = 1
eta0 = Mstar
eta1 = eta0+m1
eta2 =eta1+m2
mtilde1 = m1 * (eta0/eta1)
Mtilde1 = Mstar * (eta1/eta0)
mtilde2 = m2 * (eta1/eta2)
Mtilde2 = Mstar * (eta2/eta1)
eps = m1 * m2 / (mtilde1 + mtilde2) / Mstar
beta1 = mtilde1 / (mtilde1 + mtilde2)
beta2 = mtilde2 / (mtilde1 + mtilde2)
gamma = mtilde2/mtilde1
# Angle variable for averaging over
Q = T.dvector('Q')
# Dynamical variables:
dyvars = T.vector()
sigma1, sigma2, I1, I2, amd = [dyvars[i] for i in range(5)]
# Quadrature weights
quad_weights = T.dvector('w')
# Set lambda2=0
l2 = T.constant(0.)
l1 = -1 * k * Q
w1 = (1+s) * l2 - s * l1 - sigma1
w2 = (1+s) * l2 - s * l1 - sigma2
Gamma1 = I1
Gamma2 = I2
# Resonant semi-major axis ratio
alpha_res = ((j-k)/j)**(2/3) * ((Mstar + m1) / (Mstar+m2))**(1/3)
P0 = k * ( beta2 - beta1 * T.sqrt(alpha_res) ) / 2
P = P0 - k * (s+1/2) * amd
Ltot = beta1 * T.sqrt(alpha_res) + beta2 - amd
L1 = Ltot/2 - P / k - s * (I1 + I2)
L2 = Ltot/2 + P / k + (1 + s) * (I1 + I2)
a1 = (L1 / beta1 )**2 * eta0 / eta1
e1 = T.sqrt(1-(1-(Gamma1 / L1))**2)
a2 = (L2 / beta2 )**2 * eta1 / eta2
e2 = T.sqrt(1-(1-(Gamma2 / L2))**2)
Hkep = - eta1 * beta1 / (2 * a1) / eta0 - eta2 * beta2 / (2 * a2) / eta1
alpha = a1 / a2
ko = KeplerOp()
M1 = l1 - w1
M2 = l2 - w2
sinf1,cosf1 = ko( M1, e1 + T.zeros_like(M1) )
sinf2,cosf2 = ko( M2, e2 + T.zeros_like(M2) )
R = calc_DisturbingFunction_with_sinf_cosf(alpha,e1,e2,w1,w2,sinf1,cosf1,sinf2,cosf2)
Rav = R.dot(quad_weights)
Hpert = -eps * Rav / a2
Htot = Hkep + Hpert
######################
# Dissipative dynamics
######################
tau_alpha_0, K1, K2, p = T.dscalars(4)
sigma1dot_dis,sigma2dot_dis,I1dot_dis,I2dot_dis,amddot_dis = T.dscalars(5)
sigma1dot_dis,sigma2dot_dis = T.as_tensor(0.),T.as_tensor(0.)
# Define timescales
tau_e1 = tau_alpha_0 / K1
tau_e2 = tau_alpha_0 / K2
tau_a1_0 = -1 * tau_alpha_0 * (1+alpha_res * gamma)/ (alpha_res * gamma)
tau_a2_0 = -1 * alpha_res * gamma * tau_a1_0
tau_a1 = 1 / (1/tau_a1_0 + 2 * p * e1*e1 / tau_e1 )
tau_a2 = 1 / (1/tau_a2_0 + 2 * p * e2*e2 / tau_e2 )
# Time derivative of orbital elements
e1dot_dis = -1*e1 / tau_e1
e2dot_dis = -1*e2 / tau_e2
a1dot_dis = -1*a1 / tau_a1
a2dot_dis = -1*a2 / tau_a2
# Time derivatives of canonical variables
I1dot_dis = L1 * e1 * e1dot_dis / ( T.sqrt(1-e1*e1) ) - I1 / tau_a1 / 2
I2dot_dis = L2 * e2 * e2dot_dis / ( T.sqrt(1-e2*e2) ) - I2 / tau_a2 / 2
Pdot_dis = -1*k * ( L2 / tau_a2 - L1 / tau_a1) / 4 - k * (s + 1/2) * (I1dot_dis + I2dot_dis)
amddot_dis = Pdot_dis / T.grad(P,amd)
#####################################################
# Set parameters for compiling functions with Theano
#####################################################
# Get numerical quadrature nodes and weights
nodes,weights = np.polynomial.legendre.leggauss(N_QUAD_PTS)
# Rescale for integration interval from [-1,1] to [-pi,pi]
nodes = nodes * np.pi
weights = weights * 0.5
# 'givens' will fix some parameters of Theano functions compiled below
givens = [(Q,nodes),(quad_weights,weights)]
# 'ins' will set the inputs of Theano functions compiled below
# Note: 'extra_ins' will be passed as values of object attributes
# of the 'ResonanceEquations' class 'defined below
extra_ins = [m1,m2,j,k,tau_alpha_0,K1,K2,p]
ins = [dyvars] + extra_ins
# Define flows and jacobians.
# Conservative flow
gradHtot = T.grad(Htot,wrt=dyvars)
hessHtot = theano.gradient.hessian(Htot,wrt=dyvars)
Jtens = T.as_tensor(np.pad(getOmegaMatrix(2),(0,1),'constant'))
H_flow_vec = Jtens.dot(gradHtot)
H_flow_jac = Jtens.dot(hessHtot)
# Dissipative flow
dis_flow_vec = T.stack(sigma1dot_dis,sigma2dot_dis,I1dot_dis,I2dot_dis,amddot_dis)
dis_flow_jac = theano.gradient.jacobian(dis_flow_vec,dyvars)
# Extras
dis_timescales = [tau_a1_0,tau_a2_0,tau_e1,tau_e2]
orbels = [a1,e1,sigma1*k,a2,e2,sigma2*k]
##########################
# Compile Theano functions
##########################
if not DEBUG:
# Note that compiling can take a while
# so I've put a debugging switch here
# to skip evaluating these functions when
# desired.
Rav_fn = theano.function(
inputs=ins,
outputs=Rav,
givens=givens,
on_unused_input='ignore'
)
Hpert_av_fn = theano.function(
inputs=ins,
outputs=Hpert,
givens=givens,
on_unused_input='ignore'
)
Htot_fn = theano.function(
inputs=ins,
outputs=Htot,
givens=givens,
on_unused_input='ignore'
)
H_flow_vec_fn = theano.function(
inputs=ins,
outputs=H_flow_vec,
givens=givens,
on_unused_input='ignore'
)
H_flow_jac_fn = theano.function(
inputs=ins,
outputs=H_flow_jac,
givens=givens,
on_unused_input='ignore'
)
dis_flow_vec_fn = theano.function(
inputs=ins,
outputs=dis_flow_vec,
givens=givens,
on_unused_input='ignore'
)
dis_flow_jac_fn = theano.function(
inputs=ins,
outputs=dis_flow_jac,
givens=givens,
on_unused_input='ignore'
)
dis_timescales_fn =theano.function(
inputs=extra_ins,
outputs=dis_timescales,
givens=givens,
on_unused_input='ignore'
)
orbels_fn = theano.function(
inputs=ins,
outputs=orbels,
givens=givens,
on_unused_input='ignore'
)
else:
return [lambda x: x for _ in range(8)]
return Rav_fn,Hpert_av_fn,Htot_fn,H_flow_vec_fn,H_flow_jac_fn,dis_flow_vec_fn,dis_flow_jac_fn,dis_timescales_fn,orbels_fn
class ResonanceEquations():
"""
A class for the model describing the dynamics of a pair of planar planets
in/near a mean motion resonance.
Includes the effects of dissipation.
Attributes
----------
j : int
Together with k specifies j:j-k resonance
k : int
Order of resonance.
alpha : float
Semi-major axis ratio a_1/a_2
eps : float
Mass parameter m1*mu2 / (mu1+mu2)
m1 : float
Inner planet mass
m2 : float
Outer planet mass
"""
def __init__(self,j,k, n_quad_pts = 40, m1 = 1e-5 , m2 = 1e-5,K1=100, K2=100, tau_alpha = 1e5, p = 1):
self.j = j
self.k = k
self.m1 = m1
self.m2 = m2
self.K1 = K1
self.K2 = K2
self.tau_alpha = tau_alpha
self.p = p
self.n_quad_pts = n_quad_pts
funcs = get_compiled_theano_functions(n_quad_pts)
self._Rav_fn,self._Hpert_av_fn,self._H,self._H_flow,self._H_jac,self._dis_flow,self._dis_jac,self._times_scales,self._orbels_fn = funcs
self._omega_syn_fn, self._Hkep_fn, self._Hpert_fn,\
self._omega_syn_grad_fn, self._Hkep_grad_fn, self._Hpert_grad_fn = get_compiled_Hkep_Hpert_full()
@property
def extra_args(self):
return [self.m1,self.m2,self.j,self.k,self.tau_alpha,self.K1,self.K2,self.p]
@property
def Hpert_extra_args(self):
return [self.m1,self.m2,self.j,self.k]
@property
def mu1(self):
return self.m1 / (1 + self.m1)
@property
def mu2(self):
return self.m2 / (1 + self.m2)
@property
def beta1(self):
return self.mu1 / (self.mu1 + self.mu2)
@property
def beta1(self):
return self.mu1 / (self.mu1 + self.mu2)
@property
def beta2(self):
return self.mu2 / (self.mu1 + self.mu2)
@property
def eps(self):
return self.m1 * self.mu2 / (self.mu1 + self.mu2)
@property
def alpha(self):
return ((self.j - self.k) / self.j)**(2/3)
@property
def timescales(self):
tscales_array = self._times_scales(*self.extra_args)
return {name:scale for name,scale in zip(['tau_a1','tau_a2','tau_e1','tau_e2'],tscales_array)}
def Rav(self,z):
"""
Calculate the value of the Q-averaged disturbing function
Arguments
---------
z : ndarray
Dynamical variables
Returns
-------
float :
The value of the Hamiltonian evaluated at z.
"""
return self._Rav_fn(z,*self.extra_args)
def Hpert_osc(self,Q,z):
"""
Calculate the value of the un-averaged disturbing function
Arguments
---------
Q : float
Synodic angle
z : ndarray
Dynamical variables
Returns
-------
float :
The value of the (unaveraged) perturbation Hamiltonian evaluated at (Q,z).
"""
arg = np.concatenate(([Q],z))
Hpert_av = self._Hpert_av_fn(z,*self.extra_args)
return self._Hpert_fn(arg,*self.Hpert_extra_args) - Hpert_av
def gradHpert(self,Q,z):
"""
Calculate the value of the averaged disturbing function
Arguments
---------
Q : float
Synodic angle
z : ndarray
Dynamical variables
Returns
-------
ndarray :
The gradient of the (unaveraged) perturbation Hamiltonian evaluated at (Q,z).
"""
arg = np.concatenate(([Q],z))
return self._Hpert_grad_fn(arg,*self.Hpert_extra_args)
def H(self,z):
"""
Calculate the value of the Hamiltonian.
Arguments
---------
z : ndarray
Dynamical variables
Returns
-------
float :
The value of the Hamiltonian evaluated at z.
"""
return self._H(z,*self.extra_args)
def H_kep(self,z):
"""
Calculate the value of the Keplerian component of the Hamiltonian.
Arguments
---------
z : ndarray
Dynamical variables
Returns
-------
float :
The value of the Keplerian part of the Hamiltonian evaluated at z.
"""
Htot = self.H(z)
Hpert = self.H_pert(z)
return Htot - Hpert
def H_pert(self,z):
"""
Calculate the value of the perturbation component of the Hamiltonian.
Arguments
---------
z : ndarray
Dynamical variables
Returns
-------
float :
The value of the perturbation part of the Hamiltonian evaluated at z.
"""
Hpert = self._Hpert_av_fn(z,*self.extra_args)
return Hpert
def omega_syn(self,z):
"""
Calculate the synodic frequency dH0/dP.
Arguments
---------
z : ndarray
Dynamical variables
Returns
-------
float :
The value of dH0/dP evaluated at z.
"""
zfull = np.concatenate(([0],z))
return self._omega_syn_fn(zfull,*self.Hpert_extra_args)
def grad_omega_syn(self,z):
"""
Calculate the gradient of the synodic frequency, grad(dH0/dP).
Arguments
---------
z : ndarray
Dynamical variables
Returns
-------
ndarray :
The value of grad(dH0/dP) evaluated at z.
"""
zfull = np.concatenate(([0],z))
return self._omega_syn_grad_fn(zfull,*self.Hpert_extra_args)
def H_flow(self,z):
"""
Calculate flow induced by the Hamiltonian
.. math:
\dot{z} = \Omega \cdot \nablda_{z}H(z)
Arguments
---------
z : ndarray
Dynamical variables
Returns
-------
ndarray :
Flow vector
"""
return self._H_flow(z,*self.extra_args)
def H_flow_jac(self,z):
"""
Calculate the Jacobian of the flow induced by the Hamiltonian
.. math:
\Omega \cdot \Delta H(z)
Arguments
---------
z : ndarray
Dynamical variables
Returns
-------
ndarray :
Jacobian matrix
"""
return self._H_jac(z,*self.extra_args)
def dis_flow(self,z):
"""
Calculate flow induced by dissipative forces
.. math:
\dot{z} = f_{dis}(z)
Arguments
---------
z : ndarray
Dynamical variables
Returns
-------
ndarray :
Flow vector
"""
return self._dis_flow(z,*self.extra_args)
def dis_flow_jac(self,z):
"""
Calculate the Jacobian of the flow induced by dissipative forces
.. math:
\dot{z} = \nabla f_{dis}(z)
Arguments
---------
z : ndarray
Dynamical variables
Returns
-------
ndarray :
Flow Jacobian
"""
return self._dis_jac(z,*self.extra_args)
def flow(self,z):
"""
Calculate flow
.. math:
\dot{z} = \Omega \cdot \nablda_{z}H(z) + f_{dis}(z)
Arguments
---------
z : ndarray
Dynamical variables
Returns
-------
ndarray :
Flow vector
"""
return self._H_flow(z,*self.extra_args) + self._dis_flow(z,*self.extra_args)
def flow_jac(self,z):
"""
Calculate flow Jacobian
.. math:
\Omega \cdot \Delta H(z) + \\nabla f_{dis}(z)
Arguments
---------
z : ndarray
Dynamical variables
Returns
-------
ndarray :
Jacobian
"""
return self._H_jac(z,*self.extra_args) + self._dis_jac(z,*self.extra_args)
def dyvars_to_orbels(self,z):
r"""
Convert dynamical variables
.. math:
z = (\sigma_1,\sigma_2,I_1,I_2,{\cal C}
to orbital elements
.. math:
(a_1,e_1,\\theta_1,a_2,e_2,\\theta_2)
"""
return self._orbels_fn(z,*self.extra_args)
def orbels_to_dyvars(self,orbels):
r"""
Convert orbital elements
.. math:
(a_1,e_1,\theta_1,a_2,e_2,\theta_2)
to dynamical variables
.. math:
z = (\sigma_1,\sigma_2,I_1,I_2,{\cal C}
"""
# Total angular momentum constrained by
# Ltot = beta1 * sqrt(alpha_res) + beta2 - amd
Mstar = 1
m1 = self.m1
m2 = self.m2
k = self.k
j = self.j
alpha_res = ((j-k)/j)**(2/3) * ((Mstar + m1) / (Mstar+m2))**(1/3)
s = (j - k) / k
a1,e1,theta1,a2,e2,theta2 = orbels
L1 = self.beta1 * np.sqrt(a1)
L2 = self.beta2 * np.sqrt(a2)
Gamma1 = L1 * (1 - np.sqrt(1-e1**2))
Gamma2 = L2 * (1 - np.sqrt(1-e2**2))
I1 = Gamma1
I2 = Gamma2
P = k * (L2-L1) / 2 - k * (s+1/2)*(I1+I2)
Ltot = L1 + L2 - I1 - I2
Lcirc = self.beta1 * np.sqrt(alpha_res) + self.beta2
P0 = k * (self.beta2 - self.beta1 * np.sqrt(alpha_res)) / 2
# Now solve for re-scaling factor 'scale' and AMD value that satisfy
# scale * P = P0 - k*(s+1/2)*amd
# scale * Ltot = Lcirc - amd
# Solution is given by:
# scale = (2 k Lcirc s-k Lcirc+2 P0)/(k Ltot-2 P+2 k Ltot s)
# (2 (Lcirc P-Ltot P0))/(k Ltot-2 P+2 k Ltot s)
scale = (k * (1 + 2 * s) * Lcirc - 2 * P0) / (k * (1 + 2 * s) * Ltot - 2 * P)
amd = Lcirc - scale * Ltot
I1 *= scale
I2 *= scale
sigma1 = theta1 / k
sigma2 = theta2 / k
return np.array((sigma1,sigma2,I1,I2,amd))
def gradHkep(self,z):
zfull = np.concatenate(([0],z))
return self._Hkep_grad_fn(zfull,*self.Hpert_extra_args)
def mean_to_osculating_dyvars(self,Q,z,N = 256):
r"""
Apply perturbation theory to transfrom from the phase space coordiantes of the
averaged model to osculating phase space coordintates of the full phase space.
Assumes that the transformation is being applied at the fixed point of the
averaged model.
Arguemnts
---------
Q : float or ndarry
Value(s) of angle (lambda2-lambda1)/k at which to apply transformation.
Equilibrium points of the averaged model correspond to orbits that are
periodic in Q in the full phase space.
z : ndarray
Dynamical variables of the averaged model:
$\sigma_1,\sigma_2,I_1,I_2,AMD$
N : int, optional
Number of Q points to evaluate functions at when performing Fourier
transformation. Default is
N=256
Returns
-------
zosc : ndarray, (5,) or (M,5)
The osculating values of the dynamical varaibles for the
input Q values. The dimension of the returned variables
is set by the dimension of the input 'Q'. If Q is a
float, then z is an array of length 5. If Q is an
array of length M then zosc is an (M,5) array.
"""
omega_syn = self.omega_syn(z)
OmegaMtrx = getOmegaMatrix(2)
Omega_del2H = self.H_flow_jac(z)[:-1,:-1]
vals,S = np.linalg.eig(Omega_del2H)
S = np.transpose([S.T[i] for i in (0,2,1,3)])
s1 = (S.T @ OmegaMtrx @ S)[0,2]
s2 = (S.T @ OmegaMtrx @ S)[1,3]
S.T[0]*=1/np.sqrt(s1)
S.T[2]*=1/np.sqrt(s1)
S.T[1]*=1/np.sqrt(s2)
S.T[3]*=1/np.sqrt(s2)
Sinv = np.linalg.inv(S)
Qarr = np.atleast_1d(Q)
dchi = np.zeros((4,len(Qarr)),dtype=complex)
# Fill arrays for FFT
Qs = np.linspace(0,2*np.pi,N)
X = np.zeros((N,4),dtype=complex)
gradHkep = self.gradHkep(z)
domega_syn_dz = self.grad_omega_syn(z)[1:-1]
domega_syn_dw = S.T @ domega_syn_dz
for i,q in enumerate(Qs):
gradHosc = self.gradHpert(q,z) + gradHkep
X[i] = S.T @ (gradHosc)[1:-1]
X[i] -= self.Hpert_osc(q,z) * domega_syn_dw/omega_syn
omegas = +1*np.imag(np.diag(Sinv @ Omega_del2H @ S))[:2]
for I in range(2):
A = np.fft.fft(X[:,I])
freqs = np.fft.fftshift(np.fft.fftfreq(N)*N)
amps = np.fft.fftshift(A)/N
for l in range(1,N//2 - 1):
sig = -1j * self.k * amps[N//2+l] * np.exp(1j * freqs[N//2+l] * Qarr) / (-l*omega_syn - self.k * omegas[I])
sig +=-1j * self.k * amps[N//2-l] * np.exp(1j * freqs[N//2-l] * Qarr) / (l*omega_syn - self.k * omegas[I])
dchi[I] += sig
dchi[2] = -1j * np.conjugate(dchi[0])
dchi[3] = -1j * np.conjugate(dchi[1])
# Get AMD correction
s = (self.j - self.k) / self.k
dAMD = np.array([self.Hpert_osc(q,z) for q in Qarr]) / omega_syn / (s+1/2)
dsigmaI = np.transpose(-1 * (S @ OmegaMtrx @ dchi).T)
result = np.transpose(z + np.vstack((dsigmaI,dAMD)).T)
result = np.real(result) # trim small imaginary parts cause by numerical errors
if result.shape[1] == 1:
return result.reshape(-1)
return result
from scipy.optimize import lsq_linear
def find_equilibrium(self,guess,dissipation=False,tolerance=1e-9,max_iter=10):
"""
Use Newton's method to locate an equilibrium solution of
the equations of motion.
By default, an equilibrium of the dissipation-free equations
is sought. In this case, the AMD value of the equilibrium
solution will be equal to the value of the initially supplied
guess of dynamical variables. If dissipative terms are included
then the equilibrium will depend on the parameters K1 and K2
as well as tau_alpha.
Arguments
---------
guess : ndarray
Initial guess for the dynamical variables at the equilibrium.
dissipation : bool, optional
Whether dissipative terms are considered in the equations of
motion. Default is False.
tolerance : float, optional
Tolerance for root-finding such that solution satisfies |f(z)| < tolerance
Default value is 1E-9.
max_iter : int, optional
Maximum number of Newton's method iterations.
Default is 10.
include_dissipation : bool, optional
Include dissipative terms in the equations of motion.
Default is False
Returns
-------
zeq : ndarray
Equilibrium value of dynamical variables.
Raises
------
RuntimeError : Raises error if maximum number of Newton iterations is exceeded.
"""
if dissipation:
return self._find_dissipative_equilibrium(guess,tolerance,max_iter)
else:
return self._find_conservative_equilibrium(guess,tolerance,max_iter)
def _find_dissipative_equilibrium(self,guess,tolerance=1e-9,max_iter=10):
y = guess
lb = np.array([-np.pi,-np.pi, -1* y[2], -1 * y[3],-np.inf])
ub = np.array([np.pi,np.pi,np.inf,np.inf,np.inf])
fun = self.flow
jac = self.flow_jac
f = fun(y)
J = jac(y)
it=0
# Newton method iteration
while np.linalg.norm(f)>tolerance and it < max_iter:
# Note-- using constrained least-squares
# to avoid setting actions to negative
# values.
# The lower bounds ensure that I1 and I2 are positive quantities
lb[2:-1] = -1* y[2:-1]
dy = lsq_linear(J,-f,bounds=(lb,ub)).x
y = y + dy
f = fun(y)
J = jac(y)
it+=1
if it==max_iter:
raise RuntimeError("Max iterations reached!")
return y
def _find_conservative_equilibrium(self,guess,tolerance=1e-9,max_iter=10):
y = guess
lb = np.array([-np.pi,-np.pi, -1* y[2], -1 * y[3]])
ub = np.array([np.pi,np.pi,np.inf,np.inf])
fun = self.H_flow
jac = self.H_flow_jac
f = fun(y)[:-1]
J = jac(y)[:-1,:-1]
it=0
# Newton method iteration
while np.linalg.norm(f)>tolerance and it < max_iter:
# Note-- using constrained least-squares
# to avoid setting actions to negative
# values.
# The lower bounds ensure that I1 and I2 are positive quantities
lb[2:] = -1* y[2:-1]
dy = lsq_linear(J,-f,bounds=(lb,ub)).x
y[:-1] = y[:-1] + dy
f = fun(y)[:-1]