forked from stevertaylor/NX01
/
NX01_utils.py
1476 lines (1095 loc) · 44.5 KB
/
NX01_utils.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
"""
Created by stevertaylor
Copyright (c) 2014 Stephen R. Taylor
Code contributions by Rutger van Haasteren (piccard) and Justin Ellis (PAL/PAL2).
"""
from __future__ import division
import numpy as np
from numpy import *
import os
import math
from scipy import integrate
from scipy.integrate import odeint
from scipy import optimize
from scipy import constants as sc
from numpy import random
from scipy import special as ss
from scipy import linalg as sl
from scipy.interpolate import interp1d
from pkg_resources import resource_filename, Requirement
import numexpr as ne
import optparse
import ephem
from ephem import *
day = 24 * 3600
year = 365.25 * day
SOLAR2S = sc.G / sc.c**3 * 1.98855e30
KPC2S = sc.parsec / sc.c * 1e3
MPC2S = sc.parsec / sc.c * 1e6
def sumTermCovarianceMatrix_fast(tm, fL, gam):
"""
Calculate the power series expansion for the Hypergeometric
function in the standard power-law covariance matrix. This
version uses the Python package numexpr and is much faster
than using numpy. For now it is hardcoded to use only the
first 3 terms.
@param tm: Matrix of time lags in years
@param fL: Low frequency cutoff
@param gam: Power Law spectral index
"""
x = 2*np.pi*fL*tm
sum = ne.evaluate("1/(1-gam) - x**2/(2*(3-gam)) + x**4/(24*(5-gam))")
return sum
def makeTimeGrid(psra, psrb):
"""
Construct time-domain DM-variation
covariance matrix.
@param psra: object for pulsar 'a'
@param psrb: object for pulsar 'b'
@return: Cdm: Time-lag grid
"""
ta, tb = np.meshgrid(psra.toas, psrb.toas)
tm = np.abs(ta-tb).astype(np.float64)/365.25
return tm
def makeRedTDcov(Ared, gam_red, tm):
"""
Construct time-domain red-noise
covariance matrix.
@param Ared: Red-noise spectral amplitude
@param gam_red: Red-noise spectral slope
@param tm: time-lag matrix
@return: Cdm: Red-noise covariance matrix
"""
Tspan = tm.max()
fL = 1/(100.0*Tspan)
xgrid = 2.0*np.pi*fL*tm
Cred = ( (Ared**2.0)*(fL**(1.0-gam_red)) / (12.0*np.pi**2.0) ) * \
((ss.gamma(1.0-gam_red)*np.sin(np.pi*gam_red/2.)*ne.evaluate("xgrid**(gam_red-1.0)"))
- sumTermCovarianceMatrix_fast(tm, fL, gam_red))
Cred *= ((365.25*86400.0)**2.0)
return Cred
def makeDmTDcov(psr, Adm, gam_dm, tm):
"""
Construct time-domain DM-variation
covariance matrix.
@param psr: pulsar object
@param Adm: DM-variation spectral amplitude
@param gam_dm: DM-variation spectral slope
@param tm: time-lag matrix
@return: Cdm: DM covariance matrix
"""
Tspan = tm.max()
fL = 1/(100.0*Tspan)
xgrid = 2.0*np.pi*fL*tm
K = 2.41*10.0**(-16.0)
Dm = 1.0/(K*(psr.obs_freqs*1e6)**2.0)
DmA,DmB = np.meshgrid(Dm,Dm)
DmGrid = DmA*DmB
Cdm = ( (Adm**2.0)*(fL**(1.0-gam_dm)) / (12.0*np.pi**2.0) ) * \
((ss.gamma(1-gam_dm)*np.sin(np.pi*gam_dm/2)*ne.evaluate("xgrid**(gam_dm-1)"))
- sumTermCovarianceMatrix_fast(tm, fL, gam_dm))
Cdm *= ((365.25*86400.0)**2.0)
Cdm = np.multiply(DmGrid,Cdm)
return Cdm
def createFourierDesignmatrix_red(t, nmodes, freq=False, Tspan=None):
"""
Construct fourier design matrix from eq 11 of Lentati et al, 2013
@param t: vector of time series in seconds
@param nmodes: number of fourier coefficients to use
@param freq: option to output frequencies
@param Tspan: option to some other Tspan
@return: F: fourier design matrix
@return: f: Sampling frequencies (if freq=True)
"""
N = len(t)
F = np.zeros((N, 2*nmodes))
if Tspan is not None:
T = Tspan
else:
T = t.max() - t.min()
# define sampling frequencies
fqs = np.linspace(1/T, nmodes/T, nmodes)
# The sine/cosine modes
ct = 0
for ii in range(0, 2*nmodes-1, 2):
F[:,ii] = np.cos(2*np.pi*fqs[ct]*t)
F[:,ii+1] = np.sin(2*np.pi*fqs[ct]*t)
ct += 1
if freq:
return F, fqs
else:
return F
def createFourierDesignmatrix_dm(t, nmodes, obs_freqs, freq=False, Tspan=None):
"""
Construct fourier design matrix from eq 11 of Lentati et al, 2013
@param t: vector of time series in seconds
@param nmodes: number of fourier coefficients to use
@param freq: option to output frequencies
@param Tspan: option to some other Tspan
@param pbs_freqs: pulsar observing frequencies
@return: F: fourier design matrix
@return: f: Sampling frequencies (if freq=True)
"""
N = len(t)
F = np.zeros((N, 2*nmodes))
if Tspan is not None:
T = Tspan
else:
T = t.max() - t.min()
# define sampling frequencies
fqs = np.linspace(1/T, nmodes/T, nmodes)
# compute the DM-variation vectors
K = 2.41*10.0**(-16.0)
Dm = 1.0/(K*(obs_freqs*1e6)**2.0)
# The sine/cosine modes
ct = 0
for ii in range(0, 2*nmodes-1, 2):
F[:,ii] = np.multiply(np.cos(2*np.pi*fqs[ct]*t),Dm)
F[:,ii+1] = np.multiply(np.sin(2*np.pi*fqs[ct]*t),Dm)
ct += 1
if freq:
return F, fqs
else:
return F
def quantize_fast(times, dt=1.0, calci=False):
""" Adapted from libstempo: produce the quantisation matrix fast """
isort = np.argsort(times)
bucket_ref = [times[isort[0]]]
bucket_ind = [[isort[0]]]
for i in isort[1:]:
if times[i] - bucket_ref[-1] < dt:
bucket_ind[-1].append(i)
else:
bucket_ref.append(times[i])
bucket_ind.append([i])
t = np.array([np.mean(times[l]) for l in bucket_ind],'d')
U = np.zeros((len(times),len(bucket_ind)),'d')
for i,l in enumerate(bucket_ind):
U[l,i] = 1
rv = (t, U)
if calci:
Ui = ((1.0/np.sum(U, axis=0)) * U).T
rv = (t, U, Ui)
return rv
def quantize_split(times, flags, dt=1.0, calci=False):
"""
As quantize_fast, but now split the blocks per backend. Note: for
efficiency, this function assumes that the TOAs have been sorted by
argsortTOAs. This is _NOT_ checked.
"""
isort = np.arange(len(times))
bucket_ref = [times[isort[0]]]
bucket_flag = [flags[isort[0]]]
bucket_ind = [[isort[0]]]
for i in isort[1:]:
if times[i] - bucket_ref[-1] < dt and flags[i] == bucket_flag[-1]:
bucket_ind[-1].append(i)
else:
bucket_ref.append(times[i])
bucket_flag.append(flags[i])
bucket_ind.append([i])
t = np.array([np.mean(times[l]) for l in bucket_ind],'d')
U = np.zeros((len(times),len(bucket_ind)),'d')
for i,l in enumerate(bucket_ind):
U[l,i] = 1
rv = (t, U)
if calci:
Ui = ((1.0/np.sum(U, axis=0)) * U).T
rv = (t, U, Ui)
return rv
def argsortTOAs(toas, flags, which=None, dt=1.0):
"""
Return the sort, and the inverse sort permutations of the TOAs, for the
requested type of sorting
NOTE: This one is _not_ optimized for efficiency yet (but is done only once)
:param toas: The toas that are to be sorted
:param flags: The flags that belong to each TOA (indicates sys/backend)
:param which: Which type of sorting we will use (None, 'jitterext', 'time')
:param dt: Timescale for which to limit jitter blocks, default [10 secs]
:return: perm, perminv (sorting permutation, and inverse)
"""
if which is None:
isort = slice(None, None, None)
iisort = slice(None, None, None)
elif which == 'time':
isort = np.argsort(toas, kind='mergesort')
iisort = np.zeros(len(isort), dtype=np.int)
for ii, p in enumerate(isort):
iisort[p] = ii
elif which == 'jitterext':
tave, Umat = quantize_fast(toas, dt)
isort = np.argsort(toas, kind='mergesort')
uflagvals = list(set(flags))
for cc, col in enumerate(Umat.T):
for flagval in uflagvals:
flagmask = (flags[isort] == flagval)
if np.sum(col[isort][flagmask]) > 1:
# This observing epoch has several TOAs
colmask = col[isort].astype(np.bool)
epmsk = flagmask[colmask]
epinds = np.flatnonzero(epmsk)
if len(epinds) == epinds[-1] - epinds[0] + 1:
# Keys are exclusively in succession
pass
else:
# Sort the indices of this epoch and backend
# We need mergesort here, because it is stable
# (A stable sort keeps items with the same key in the
# same relative order. )
episort = np.argsort(flagmask[colmask], kind='mergesort')
isort[colmask] = isort[colmask][episort]
else:
# Only one element, always ok
pass
# Now that we have a correct permutation, also construct the inverse
iisort = np.zeros(len(isort), dtype=np.int)
for ii, p in enumerate(isort):
iisort[p] = ii
else:
isort, iisort = np.arange(len(toas)), np.arange(len(toas))
return isort, iisort
def checkTOAsort(toas, flags, which=None, dt=1.0):
"""
Check whether the TOAs are indeed sorted as they should be according to the
definition in argsortTOAs
:param toas: The toas that are supposed to be already sorted
:param flags: The flags that belong to each TOA (indicates sys/backend)
:param which: Which type of sorting we will check (None, 'jitterext', 'time')
:param dt: Timescale for which to limit jitter blocks, default [10 secs]
:return: True/False
"""
rv = True
if which is None:
isort = slice(None, None, None)
iisort = slice(None, None, None)
elif which == 'time':
isort = np.argsort(toas, kind='mergesort')
if not np.all(isort == np.arange(len(isort))):
rv = False
elif which == 'jitterext':
tave, Umat = quantize_fast(toas, dt)
#isort = np.argsort(toas, kind='mergesort')
isort = np.arange(len(toas))
uflagvals = list(set(flags))
for cc, col in enumerate(Umat.T):
for flagval in uflagvals:
flagmask = (flags[isort] == flagval)
if np.sum(col[isort][flagmask]) > 1:
# This observing epoch has several TOAs
colmask = col[isort].astype(np.bool)
epmsk = flagmask[colmask]
epinds = np.flatnonzero(epmsk)
if len(epinds) == epinds[-1] - epinds[0] + 1:
# Keys are exclusively in succession
pass
else:
# Keys are not sorted for this epoch/flag
rv = False
else:
# Only one element, always ok
pass
else:
pass
return rv
def checkquant(U, flags, uflagvals=None):
"""
Check the quantization matrix for consistency with the flags
:param U: quantization matrix
:param flags: the flags of the TOAs
:param uflagvals: subset of flags that are not ignored
:return: True/False, whether or not consistent
The quantization matrix is checked for three kinds of consistency:
- Every quantization epoch has more than one observation
- No quantization epoch has no observations
- Only one flag is allowed per epoch
"""
if uflagvals is None:
uflagvals = list(set(flags))
rv = True
collisioncheck = np.zeros((U.shape[1], len(uflagvals)), dtype=np.int)
for ii, flagval in enumerate(uflagvals):
flagmask = (flags == flagval)
Umat = U[flagmask, :]
simepoch = np.sum(Umat, axis=0)
if np.all(simepoch <= 1) and not np.all(simepoch == 0):
rv = False
#raise ValueError("quantization matrix contains non-jitter-style data")
collisioncheck[:, ii] = simepoch
# Check continuity of the columns
for cc, col in enumerate(Umat.T):
if np.sum(col > 2):
# More than one TOA for this flag/epoch
epinds = np.flatnonzero(col)
if len(epinds) != epinds[-1] - epinds[0] + 1:
rv = False
print("WARNING: checkquant found non-continuous blocks")
#raise ValueError("quantization matrix epochs not continuous")
epochflags = np.sum(collisioncheck > 0, axis=1)
if np.any(epochflags > 1):
rv = False
print("WARNING: checkquant found multiple backends for an epoch")
print epochflags
#raise ValueError("Some observing epochs include multiple backends")
if np.any(epochflags < 1):
rv = False
print("WARNING: checkquant found epochs without observations (eflags)")
#raise ValueError("Some observing epochs include no observations... ???")
obsum = np.sum(U, axis=0)
if np.any(obsum < 1):
rv = False
print("WARNING: checkquant found epochs without observations (all)")
#raise ValueError("Some observing epochs include no observations... ???")
return rv
def quant2ind(U):
"""
Convert the quantization matrix to an indices matrix for fast use in the
jitter likelihoods
:param U: quantization matrix
:return: Index (basic slicing) version of the quantization matrix
This function assumes that the TOAs have been properly sorted according to
the proper function argsortTOAs above. Checks on the continuity of U are not
performed
"""
inds = np.zeros((U.shape[1], 2), dtype=np.int)
for cc, col in enumerate(U.T):
epinds = np.flatnonzero(col)
inds[cc, 0] = epinds[0]
inds[cc, 1] = epinds[-1]+1
return inds
def quantreduce(U, eat, flags, calci=False):
"""
Reduce the quantization matrix by removing the observing epochs that do not
require any jitter parameters.
:param U: quantization matrix
:param eat: Epoch-averaged toas
:param flags: the flags of the TOAs
:param calci: Calculate pseudo-inverse yes/no
:return newU, jflags (flags that need jitter)
"""
uflagvals = list(set(flags))
incepoch = np.zeros(U.shape[1], dtype=np.bool)
jflags = []
for ii, flagval in enumerate(uflagvals):
flagmask = (flags == flagval)
Umat = U[flagmask, :]
ecnt = np.sum(Umat, axis=0)
incepoch = np.logical_or(incepoch, ecnt>1)
if np.any(ecnt > 1):
jflags.append(flagval)
Un = U[:, incepoch]
eatn = eat[incepoch]
if calci:
Ui = ((1.0/np.sum(Un, axis=0)) * Un).T
rv = (Un, Ui, eatn, jflags)
else:
rv = (Un, eatn, jflags)
return rv
def dailyAve(times, res, err, ecorr, dt=1, flags=None):
# Does not work yet in NX01"
isort = np.argsort(times)
bucket_ref = [times[isort[0]]]
bucket_ind = [[isort[0]]]
for i in isort[1:]:
if times[i] - bucket_ref[-1] < dt:
bucket_ind[-1].append(i)
else:
bucket_ref.append(times[i])
bucket_ind.append([i])
avetoas = np.array([np.mean(times[l]) for l in bucket_ind],'d')
if flags is not None:
aveflags = np.array([flags[l[0]] for l in bucket_ind])
aveerr = np.zeros(len(bucket_ind))
averes = np.zeros(len(bucket_ind))
for i,l in enumerate(bucket_ind):
M = np.ones(len(l))
C = np.diag(err[l]**2) + np.ones((len(l), len(l))) * ecorr[l[0]]
avr = 1/np.dot(M, np.dot(np.linalg.inv(C), M))
aveerr[i] = np.sqrt(avr)
averes[i] = avr * np.dot(M, np.dot(np.linalg.inv(C), res[l]))
if flags is not None:
return avetoas, averes, aveerr, aveflags
else:
return avetoas, aveerr, averes
def make_ecc_interpolant():
"""
Make interpolation function from eccentricity file to
determine number of harmonics to use for a given
eccentricity.
:returns: interpolant
"""
pth = resource_filename(Requirement.parse('libstempo'),
'libstempo/ecc_vs_nharm.txt')
fil = np.loadtxt(pth)
return interp1d(fil[:,0], fil[:,1])
# get interpolant for eccentric binaries
ecc_interp = make_ecc_interpolant()
def get_edot(F, mc, e):
"""
Compute eccentricity derivative from Taylor et al. (2015)
:param F: Orbital frequency [Hz]
:param mc: Chirp mass of binary [Solar Mass]
:param e: Eccentricity of binary
:returns: de/dt
"""
# chirp mass
mc *= SOLAR2S
dedt = -304/(15*mc) * (2*np.pi*mc*F)**(8/3) * e * \
(1 + 121/304*e**2) / ((1-e**2)**(5/2))
return dedt
def get_Fdot(F, mc, e):
"""
Compute frequency derivative from Taylor et al. (2015)
:param F: Orbital frequency [Hz]
:param mc: Chirp mass of binary [Solar Mass]
:param e: Eccentricity of binary
:returns: dF/dt
"""
# chirp mass
mc *= SOLAR2S
dFdt = 48 / (5*np.pi*mc**2) * (2*np.pi*mc*F)**(11/3) * \
(1 + 73/24*e**2 + 37/96*e**4) / ((1-e**2)**(7/2))
return dFdt
def get_gammadot(F, mc, q, e):
"""
Compute gamma dot from Barack and Cutler (2004)
:param F: Orbital frequency [Hz]
:param mc: Chirp mass of binary [Solar Mass]
:param q: Mass ratio of binary
:param e: Eccentricity of binary
:returns: dgamma/dt
"""
# chirp mass
mc *= SOLAR2S
#total mass
m = (((1+q)**2)/q)**(3/5) * mc
dgdt = 6*np.pi*F * (2*np.pi*F*m)**(2/3) / (1-e**2) * \
(1 + 0.25*(2*np.pi*F*m)**(2/3)/(1-e**2)*(26-15*e**2))
return dgdt
def get_coupled_ecc_eqns(y, t, mc, q):
"""
Computes the coupled system of differential
equations from Peters (1964) and Barack &
Cutler (2004). This is a system of three variables:
F: Orbital frequency [Hz]
e: Orbital eccentricity
gamma: Angle of precession of periastron [rad]
phase0: Orbital phase [rad]
:param y: Vector of input parameters [F, e, gamma]
:param t: Time [s]
:param mc: Chirp mass of binary [Solar Mass]
:param q: Mass ratio of binary
:returns: array of derivatives [dF/dt, de/dt, dgamma/dt, dphase/dt]
"""
F = y[0]
e = y[1]
gamma = y[2]
phase = y[3]
#total mass
m = (((1+q)**2)/q)**(3/5) * mc
dFdt = get_Fdot(F, mc, e)
dedt = get_edot(F, mc, e)
dgdt = get_gammadot(F, mc, q, e)
dphasedt = 2*np.pi*F
return np.array([dFdt, dedt, dgdt, dphasedt])
def solve_coupled_ecc_solution(F0, e0, gamma0, phase0, mc, q, t):
"""
Compute the solution to the coupled system of equations
from from Peters (1964) and Barack & Cutler (2004) at
a given time.
:param F0: Initial orbital frequency [Hz]
:param e0: Initial orbital eccentricity
:param gamma0: Initial angle of precession of periastron [rad]
:param mc: Chirp mass of binary [Solar Mass]
:param q: Mass ratio of binary
:param t: Time at which to evaluate solution [s]
:returns: (F(t), e(t), gamma(t), phase(t))
"""
y0 = np.array([F0, e0, gamma0, phase0])
y, infodict = odeint(get_coupled_ecc_eqns, y0, t, args=(mc,q), full_output=True)
if infodict['message'] == 'Integration successful.':
ret = y
else:
ret = 0
return ret
def get_an(n, mc, dl, F, e):
"""
Compute a_n from Eq. 22 of Taylor et al. (2015).
:param n: Harmonic number
:param mc: Chirp mass of binary [Solar Mass]
:param dl: Luminosity distance [Mpc]
:param F: Orbital frequency of binary [Hz]
:param e: Orbital Eccentricity
:returns: a_n
"""
# convert to seconds
mc *= SOLAR2S
dl *= MPC2S
omega = 2 * np.pi * F
amp = n * mc**(5/3) * omega**(2/3) / dl
ret = -amp * (ss.jn(n-2,n*e) - 2*e*ss.jn(n-1,n*e) +
(2/n)*ss.jn(n,n*e) + 2*e*ss.jn(n+1,n*e) -
ss.jn(n+2,n*e))
return ret
def get_bn(n, mc, dl, F, e):
"""
Compute b_n from Eq. 22 of Taylor et al. (2015).
:param n: Harmonic number
:param mc: Chirp mass of binary [Solar Mass]
:param dl: Luminosity distance [Mpc]
:param F: Orbital frequency of binary [Hz]
:param e: Orbital Eccentricity
:returns: b_n
"""
# convert to seconds
mc *= SOLAR2S
dl *= MPC2S
omega = 2 * np.pi * F
amp = n * mc**(5/3) * omega**(2/3) / dl
ret = -amp * np.sqrt(1-e**2) *(ss.jn(n-2,n*e) - 2*ss.jn(n,n*e) +
ss.jn(n+2,n*e))
return ret
def get_cn(n, mc, dl, F, e):
"""
Compute c_n from Eq. 22 of Taylor et al. (2015).
:param n: Harmonic number
:param mc: Chirp mass of binary [Solar Mass]
:param dl: Luminosity distance [Mpc]
:param F: Orbital frequency of binary [Hz]
:param e: Orbital Eccentricity
:returns: c_n
"""
# convert to seconds
mc *= SOLAR2S
dl *= MPC2S
omega = 2 * np.pi * F
amp = 2 * mc**(5/3) * omega**(2/3) / dl
ret = amp * ss.jn(n,n*e) / (n * omega)
return ret
def calculate_splus_scross(nmax, mc, dl, F, e, t, l0, gamma, gammadot, inc):
"""
Calculate splus and scross summed over all harmonics.
This waveform differs slightly from that in Taylor et al (2015)
in that it includes the time dependence of the advance of periastron.
:param nmax: Total number of harmonics to use
:param mc: Chirp mass of binary [Solar Mass]
:param dl: Luminosity distance [Mpc]
:param F: Orbital frequency of binary [Hz]
:param e: Orbital Eccentricity
:param t: TOAs [s]
:param l0: Initial eccentric anomoly [rad]
:param gamma: Angle of periastron advance [rad]
:param gammadot: Time derivative of angle of periastron advance [rad/s]
:param inc: Inclination angle [rad]
"""
n = np.arange(1, nmax)
# time dependent amplitudes
an = get_an(n, mc, dl, F, e)
bn = get_bn(n, mc, dl, F, e)
cn = get_cn(n, mc, dl, F, e)
# time dependent terms
omega = 2*np.pi*F
gt = gamma + gammadot * t
lt = l0 + omega * t
# tiled phase
phase1 = n * np.tile(lt, (nmax-1,1)).T
phase2 = np.tile(gt, (nmax-1,1)).T
#phasep = phase1 + 2*phase2
#phasem = phase1 - 2*phase2
# intermediate terms
#sp = np.sin(phasem)/(n*omega-2*gammadot) + \
# np.sin(phasep)/(n*omega+2*gammadot)
#sm = np.sin(phasem)/(n*omega-2*gammadot) - \
# np.sin(phasep)/(n*omega+2*gammadot)
#cp = np.cos(phasem)/(n*omega-2*gammadot) + \
# np.cos(phasep)/(n*omega+2*gammadot)
#cm = np.cos(phasem)/(n*omega-2*gammadot) - \
# np.cos(phasep)/(n*omega+2*gammadot)
#
#
#splus_n = -0.5 * (1+np.cos(inc)**2) * (an*sp - bn*sm) + \
# (1-np.cos(inc)**2)*cn * np.sin(phase1)
#scross_n = np.cos(inc) * (an*cm - bn*cp)
sinp1 = np.sin(phase1)
cosp1 = np.cos(phase1)
sinp2 = np.sin(2*phase2)
cosp2 = np.cos(2*phase2)
sinpp = sinp1*cosp2 + cosp1*sinp2
cospp = cosp1*cosp2 - sinp1*sinp2
sinpm = sinp1*cosp2 - cosp1*sinp2
cospm = cosp1*cosp2 + sinp1*sinp2
# intermediate terms
sp = sinpm/(n*omega-2*gammadot) + \
sinpp/(n*omega+2*gammadot)
sm = sinpm/(n*omega-2*gammadot) - \
sinpp/(n*omega+2*gammadot)
cp = cospm/(n*omega-2*gammadot) + \
cospp/(n*omega+2*gammadot)
cm = cospm/(n*omega-2*gammadot) - \
cospp/(n*omega+2*gammadot)
splus_n = -0.5 * (1+np.cos(inc)**2) * (an*sp - bn*sm) + \
(1-np.cos(inc)**2)*cn * sinp1
scross_n = np.cos(inc) * (an*cm - bn*cp)
return np.sum(splus_n, axis=1), np.sum(scross_n, axis=1)
def fplus_fcross(psr, gwtheta, gwphi):
"""
Compute gravitational-wave quadrupolar antenna pattern.
:param psr: pulsar object
:param gwtheta: Polar angle of GW source in celestial coords [radians]
:param gwphi: Azimuthal angle of GW source in celestial coords [radians]
:returns: fplus, fcross
"""
# define variable for later use
cosgwtheta, cosgwphi = np.cos(gwtheta), np.cos(gwphi)
singwtheta, singwphi = np.sin(gwtheta), np.sin(gwphi)
# unit vectors to GW source
m = np.array([singwphi, -cosgwphi, 0.0])
n = np.array([-cosgwtheta*cosgwphi, -cosgwtheta*singwphi, singwtheta])
omhat = np.array([-singwtheta*cosgwphi, -singwtheta*singwphi, -cosgwtheta])
# pulsar location
ptheta = np.pi/2 - psr.psr_locs[1]
pphi = psr.psr_locs[0]
# use definition from Sesana et al 2010 and Ellis et al 2012
phat = np.array([np.sin(ptheta)*np.cos(pphi), np.sin(ptheta)*np.sin(pphi),\
np.cos(ptheta)])
fplus = 0.5 * (np.dot(m, phat)**2 - np.dot(n, phat)**2) / (1+np.dot(omhat, phat))
fcross = (np.dot(m, phat)*np.dot(n, phat)) / (1 + np.dot(omhat, phat))
return fplus, fcross
def ecc_cgw_signal(psr, gwtheta, gwphi, mc, dist, F, inc, psi, gamma0,
e0, l0, q, nmax=100, nset=None, pd=None, periEv=True,
psrTerm=False, tref=0, check=False, useFile=True,
epochTOAs=False, dummy_toas=None):
"""
Simulate GW from eccentric SMBHB. Waveform models from
Taylor et al. (2015) and Barack and Cutler (2004).
WARNING: This residual waveform is only accurate if the
GW frequency is not significantly evolving over the
observation time of the pulsar.
:param psr: pulsar object
:param gwtheta: Polar angle of GW source in celestial coords [radians]
:param gwphi: Azimuthal angle of GW source in celestial coords [radians]
:param mc: Chirp mass of SMBMB [solar masses]
:param dist: Luminosity distance to SMBMB [Mpc]
:param F: Orbital frequency of SMBHB [Hz]
:param inc: Inclination of GW source [radians]
:param psi: Polarization of GW source [radians]
:param gamma0: Initial angle of periastron [radians]
:param e0: Initial eccentricity of SMBHB
:param l0: Initial mean anomaly [radians]
:param q: Mass ratio of SMBHB
:param nmax: Number of harmonics to use in waveform decomposition
:param nset: Fix the number of harmonics to be injected
:param pd: Pulsar distance [kpc]
:param periEv: Evolve the position of periapsis [boolean]
:param psrTerm: Option to include pulsar term [boolean]
:param tref: Fiducial time at which initial parameters are referenced [s]
:param check: Check if frequency evolves significantly over obs. time
:param useFile: Use pre-computed table of number of harmonics vs eccentricity
:returns: Vector of induced residuals
"""
# define variable for later use
cosgwtheta, cosgwphi = np.cos(gwtheta), np.cos(gwphi)
singwtheta, singwphi = np.sin(gwtheta), np.sin(gwphi)
sin2psi, cos2psi = np.sin(2*psi), np.cos(2*psi)
# unit vectors to GW source
m = np.array([singwphi, -cosgwphi, 0.0])
n = np.array([-cosgwtheta*cosgwphi, -cosgwtheta*singwphi, singwtheta])
omhat = np.array([-singwtheta*cosgwphi, -singwtheta*singwphi, -cosgwtheta])
# pulsar location
ptheta = np.pi/2 - psr.psr_locs[1]
pphi = psr.psr_locs[0]
# use definition from Sesana et al 2010 and Ellis et al 2012
phat = np.array([np.sin(ptheta)*np.cos(pphi), np.sin(ptheta)*np.sin(pphi),\
np.cos(ptheta)])
fplus = 0.5 * (np.dot(m, phat)**2 - np.dot(n, phat)**2) / (1+np.dot(omhat, phat))
fcross = (np.dot(m, phat)*np.dot(n, phat)) / (1 + np.dot(omhat, phat))
cosMu = -np.dot(omhat, phat)
# get values from pulsar object
if dummy_toas is None:
if epochTOAs:
toas = (psr.detsig_avetoas.copy() - tref)*86400.0
elif not epochTOAs:
toas = (psr.toas.copy() - tref)*86400.0
elif dummy_toas is not None:
toas = (dummy_toas.copy() - tref)*86400.0
if check:
# check that frequency is not evolving significantly over obs. time
y = solve_coupled_ecc_solution(F, e0, gamma0, l0, mc, q,
np.array([0.0,toas.max()]))
# initial and final values over observation time
Fc0, ec0, gc0, phic0 = y[0,:]
Fc1, ec1, gc1, phic1 = y[-1,:]
# observation time
Tobs = 1/(toas.max()-toas.min())
if np.abs(Fc0-Fc1) > 1/Tobs:
print('WARNING: Frequency is evolving over more than one frequency bin.')
print('F0 = {0}, F1 = {1}, delta f = {2}'.format(Fc0, Fc1, 1/Tobs))
# get gammadot for earth term
if not periEv:
gammadot = 0.0
else:
gammadot = get_gammadot(F, mc, q, e0)
if nset is not None:
nharm = nset
elif useFile:
if e0 > 0.001 and e0 < 0.999:
nharm = min(int(ecc_interp(e0)), nmax) + 1
elif e0 <= 0.001:
nharm = 3