-
Notifications
You must be signed in to change notification settings - Fork 490
/
obsdata.py
3622 lines (2854 loc) · 144 KB
/
obsdata.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
# obsdata.py
# a interferometric observation class
#
# Copyright (C) 2018 Andrew Chael
#
# This program is free software: you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# This program is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU General Public License for more details.
#
# You should have received a copy of the GNU General Public License
# along with this program. If not, see <http://www.gnu.org/licenses/>.
from __future__ import division
from __future__ import print_function
from builtins import str
from builtins import range
from builtins import object
import string, copy
import numpy as np
import numpy.lib.recfunctions as rec
import matplotlib as mpl
import matplotlib.pyplot as plt
import scipy.optimize as opt
import itertools as it
import sys
import copy
import ehtim.image
import ehtim.io.save
import ehtim.io.load
import ehtim.observing.obs_simulate as simobs
from ehtim.const_def import *
from ehtim.observing.obs_helpers import *
from ehtim.statistics.dataframes import *
from ehtim.statistics.stats import *
import scipy.spatial as spatial
import warnings
warnings.filterwarnings("ignore", message="Casting complex values to real discards the imaginary part")
RAPOS = 0
DECPOS = 1
RFPOS = 2
BWPOS = 3
DATPOS = 4
TARRPOS = 5
##################################################################################################
# Obsdata object
##################################################################################################
class Obsdata(object):
"""A polarimetric VLBI observation of visibility amplitudes and phases (in Jy).
Attributes:
source (str): The source name
ra (float): The source Right Ascension in fractional hours
dec (float): The source declination in fractional degrees
mjd (int): The integer MJD of the observation
tstart (float): The start time of the observation in hours
tstop (float): The end time of the observation in hours
rf (float): The observation frequency in Hz
bw (float): The observation bandwidth in Hz
timetype (str): How to interpret tstart and tstop; either 'GMST' or 'UTC'
polrep (str): polarization representation, either 'stokes' or 'circ'
tarr (numpy.recarray): The array of telescope data with datatype DTARR
tkey (dict): A dictionary of rows in the tarr for each site name
data (numpy.recarray): the basic data with datatype DTPOL_STOKES or DTPOL_CIRC
scantable (numpy.recarray): The array of scan information
ampcal (bool): True if amplitudes calibrated
phasecal (bool): True if phases calibrated
opacitycal (bool): True if time-dependent opacities correctly accounted for in sigmas
frcal (bool): True if feed rotation calibrated out of visibilities
dcal (bool): True if D terms calibrated out of visibilities
amp (numpy.recarray): An array of saved (averaged) visibility amplitudes
bispec (numpy.recarray): An array of saved (averaged) bispectra
cphase (numpy.recarray): An array of saved (averaged) closure phases
camp (numpy.recarray): An array of saved (averaged) closure amplitudes
logcamp (numpy.recarray): An array of saved (averaged) log closure amplitudes
"""
def __init__(self, ra, dec, rf, bw, datatable, tarr, scantable=None,
polrep='stokes', source=SOURCE_DEFAULT, mjd=MJD_DEFAULT, timetype='UTC',
ampcal=True, phasecal=True, opacitycal=True, dcal=True, frcal=True):
"""A polarimetric VLBI observation of visibility amplitudes and phases (in Jy).
Args:
ra (float): The source Right Ascension in fractional hours
dec (float): The source declination in fractional degrees
rf (float): The observation frequency in Hz
bw (float): The observation bandwidth in Hz
datatable (numpy.recarray): the basic data with datatype DTPOL_STOKES or DTPOL_CIRC
tarr (numpy.recarray): The array of telescope data with datatype DTARR
scantable (numpy.recarray): The array of scan information
polrep (str): polarization representation, either 'stokes' or 'circ'
source (str): The source name
mjd (int): The integer MJD of the observation
timetype (str): How to interpret tstart and tstop; either 'GMST' or 'UTC'
ampcal (bool): True if amplitudes calibrated
phasecal (bool): True if phases calibrated
opacitycal (bool): True if time-dependent opacities correctly accounted for in sigmas
frcal (bool): True if feed rotation calibrated out of visibilities
dcal (bool): True if D terms calibrated out of visibilities
Returns:
obsdata (Obsdata): an Obsdata object
"""
if len(datatable) == 0:
raise Exception("No data in input table!")
if not (datatable.dtype in [DTPOL_STOKES, DTPOL_CIRC]):
raise Exception("Data table dtype should be DTPOL_STOKES or DTPOL_CIRC")
# Polarization Representation
if polrep=='stokes':
self.polrep = 'stokes'
self.poldict = POLDICT_STOKES
self.poltype = DTPOL_STOKES
elif polrep=='circ':
self.polrep = 'circ'
self.poldict = POLDICT_CIRC
self.poltype = DTPOL_CIRC
else:
raise Exception("only 'stokes' and 'circ' are supported polreps!")
# Set the various observation parameters
self.source = str(source)
self.ra = float(ra)
self.dec = float(dec)
self.rf = float(rf)
self.bw = float(bw)
self.ampcal = bool(ampcal)
self.phasecal = bool(phasecal)
self.opacitycal = bool(opacitycal)
self.dcal = bool(dcal)
self.frcal = bool(frcal)
if timetype not in ['GMST', 'UTC']:
raise Exception("timetype must by 'GMST' or 'UTC'")
self.timetype = timetype
# Save the data
self.data = datatable
self.scans = scantable
# Telescope array: default ordering is by sefd
self.tarr = tarr
self.tkey = {self.tarr[i]['site']: i for i in range(len(self.tarr))}
if np.any(self.tarr['sefdr']!=0) or np.any(self.tarr['sefdl']!=0):
self.reorder_tarr_sefd()
self.reorder_baselines()
# Get tstart, mjd and tstop
times = self.unpack(['time'])['time']
self.tstart = times[0]
self.mjd = int(mjd)
self.tstop = times[-1]
if self.tstop < self.tstart:
self.tstop += 24.0
# Saved closure quantity arrays
# TODO should we precompute these?
self.amp=None
self.bispec=None
self.cphase=None
self.camp=None
self.logcamp=None
def obsdata_args(self):
""""Copy arguments for making a new obsdata argument into a list and dictonary
"""
arglist = [self.ra, self.dec, self.rf, self.bw, self.data, self.tarr]
argdict = {'scantable': self.scans, 'polrep': self.polrep, 'source':self.source,
'mjd': self.mjd, 'timetype': self.timetype,
'ampcal':self.ampcal, 'phasecal':self.phasecal, 'opacitycal':self.opacitycal,
'dcal':self.dcal, 'frcal':self.frcal}
return (arglist, argdict)
def copy(self):
"""Copy the observation object.
Args:
Returns:
(Obsdata): a copy of the Obsdata object.
"""
newobs = copy.deepcopy(self)
#arglist, argdict = self.obsdata_args()
#newobs = Obsdata(*arglist, **argdict)
## copy over any precomputed tables
#newobs.amp = self.amp
#newobs.bispec = self.bispec
#newobs.cphase = self.cphase
#newobs.camp = self.camp
#newobs.logcamp = self.logcamp
return newobs
def switch_polrep(self, polrep_out='stokes', allow_singlepol=True, singlepol_hand='R'):
"""Return a new observation with the polarization representation changed
Args:
polrep_out (str): the polrep of the output data
allow_singlepol (bool): If True, treat single-polarization data as Stokes I when converting from 'circ' polrep to 'stokes'
singlepol_hand (str): 'R' or 'L'; determines which parallel-hand is assumed when converting 'stokes' to 'circ' when only I is present
Returns:
(Obsdata): new Obsdata object with potentially different polrep
"""
if polrep_out not in ['stokes','circ']:
raise Exception("polrep_out must be either 'stokes' or 'circ'")
if polrep_out==self.polrep:
return self.copy()
elif polrep_out=='stokes': #circ -> stokes
data = np.empty(len(self.data), dtype=DTPOL_STOKES)
rrmask = np.isnan(self.data['rrvis'])
llmask = np.isnan(self.data['llvis'])
for f in DTPOL_STOKES:
f = f[0]
if f in ['time','tint','t1', 't2', 'tau1', 'tau2','u','v']:
data[f] = self.data[f]
elif f=='vis':
data[f] = 0.5*(self.data['rrvis'] + self.data['llvis'])
elif f=='qvis':
data[f] = 0.5*(self.data['lrvis'] + self.data['rlvis'])
elif f=='uvis':
data[f] = 0.5j*(self.data['lrvis'] - self.data['rlvis'])
elif f=='vvis':
data[f] = 0.5*(self.data['rrvis'] - self.data['llvis'])
elif f in ['sigma','vsigma']:
data[f] = 0.5*np.sqrt(self.data['rrsigma']**2 + self.data['llsigma']**2)
elif f in ['qsigma','usigma']:
data[f] = 0.5*np.sqrt(self.data['rlsigma']**2 + self.data['lrsigma']**2)
if allow_singlepol:
# In cases where only one polarization is present, use it as an estimator for Stokes I
data['vis'][rrmask] = self.data['llvis'][rrmask]
data['sigma'][rrmask] = self.data['llsigma'][rrmask]
data['vis'][llmask] = self.data['rrvis'][llmask]
data['sigma'][llmask] = self.data['rrsigma'][llmask]
elif polrep_out=='circ': #stokes -> circ
data = np.empty(len(self.data), dtype=DTPOL_CIRC)
Vmask = np.isnan(self.data['vvis'])
for f in DTPOL_CIRC:
f = f[0]
if f in ['time','tint','t1', 't2', 'tau1', 'tau2','u','v']:
data[f] = self.data[f]
elif f=='rrvis':
data[f] = (self.data['vis'] + self.data['vvis'])
elif f=='llvis':
data[f] = (self.data['vis'] - self.data['vvis'])
elif f=='rlvis':
data[f] = (self.data['qvis'] + 1j*self.data['uvis'])
elif f=='lrvis':
data[f] = (self.data['qvis'] - 1j*self.data['uvis'])
elif f in ['rrsigma','llsigma']:
data[f] = np.sqrt(self.data['sigma']**2 + self.data['vsigma']**2)
elif f in ['rlsigma','lrsigma']:
data[f] = np.sqrt(self.data['qsigma']**2 + self.data['usigma']**2)
if allow_singlepol:
# In cases where only Stokes I is present, copy it to a specified parallel-hand
prefix = singlepol_hand.lower() + singlepol_hand.lower() # rr or ll
if prefix not in ['rr','ll']:
raise Exception('singlepol_hand must be R or L')
data[prefix + 'vis'][Vmask] = self.data['vis'][Vmask]
data[prefix + 'sigma'][Vmask] = self.data['sigma'][Vmask]
arglist, argdict = self.obsdata_args()
arglist[DATPOS] = data
argdict['polrep'] = polrep_out
newobs = Obsdata(*arglist, **argdict)
return newobs
def reorder_baselines(self):
"""Reorder baselines to match uvfits convention, based on the telescope array ordering
"""
# Time partition the datatable
datatable = self.data.copy()
datalist = []
for key, group in it.groupby(datatable, lambda x: x['time']):
datalist.append(np.array([obs for obs in group]))
# Remove conjugate baselines
obsdata = []
for tlist in datalist:
blpairs = []
for dat in tlist:
if not (set((dat['t1'], dat['t2']))) in blpairs:
# Reverse the baseline in the right order for uvfits:
if(self.tkey[dat['t2']] < self.tkey[dat['t1']]):
(dat['t1'], dat['t2']) = (dat['t2'], dat['t1'])
(dat['tau1'], dat['tau2']) = (dat['tau2'], dat['tau1'])
dat['u'] = -dat['u']
dat['v'] = -dat['v']
if self.polrep=='stokes':
dat['vis'] = np.conj(dat['vis'])
dat['qvis'] = np.conj(dat['qvis'])
dat['uvis'] = np.conj(dat['uvis'])
dat['vvis'] = np.conj(dat['vvis'])
elif self.polrep=='circ':
dat['rrvis'] = np.conj(dat['rrvis'])
dat['llvis'] = np.conj(dat['llvis'])
#must switch l & r !!
rl = dat['rlvis'].copy()
lr = dat['lrvis'].copy()
dat['rlvis'] = np.conj(lr)
dat['lrvis'] = np.conj(rl)
else:
raise Exception("polrep must be either 'stokes' or 'circ'")
# Append the data point
blpairs.append(set((dat['t1'],dat['t2'])))
obsdata.append(dat)
obsdata = np.array(obsdata, dtype=self.poltype)
# Timesort data
obsdata = obsdata[np.argsort(obsdata, order=['time','t1'])]
# Save the data
self.data = obsdata
return
def reorder_tarr_sefd(self):
"""Reorder the telescope array by SEFD minimal to maximum.
"""
sorted_list = sorted(self.tarr, key=lambda x: np.sqrt(x['sefdr']**2 + x['sefdl']**2))
self.tarr = np.array(sorted_list,dtype=DTARR)
self.tkey = {self.tarr[i]['site']: i for i in range(len(self.tarr))}
self.reorder_baselines()
return
def reorder_tarr_snr(self):
"""Reorder the telescope array by median SNR maximal to minimal.
"""
snr = self.unpack(['t1','t2','snr'])
snr_median = [np.median(snr[(snr['t1']==site) + (snr['t2']==site)]['snr']) for site in self.tarr['site']]
idx = np.argsort(snr_median)[::-1]
self.tarr = self.tarr[idx]
self.tkey = {self.tarr[i]['site']: i for i in range(len(self.tarr))}
self.reorder_baselines()
return
def reorder_tarr_random(self):
"""Randomly reorder the telescope array.
"""
idx = np.arange(len(self.tarr))
np.random.shuffle(idx)
self.tarr = self.tarr[idx]
self.tkey = {self.tarr[i]['site']: i for i in range(len(self.tarr))}
self.reorder_baselines()
return
def data_conj(self):
"""Make a data array including all conjugate baselines.
Args:
Returns:
(numpy.recarray): a copy of the Obsdata.data table including all conjugate baselines.
"""
data = np.empty(2*len(self.data), dtype=self.poltype)
# Add the conjugate baseline data
for f in self.poltype:
f = f[0]
if f in ['t1', 't2', 'tau1', 'tau2']:
if f[-1]=='1': f2 = f[:-1]+'2'
else: f2 = f[:-1]+'1'
data[f] = np.hstack((self.data[f], self.data[f2]))
elif f in ['u','v']:
data[f] = np.hstack((self.data[f], -self.data[f]))
elif f in [self.poldict['vis1'],self.poldict['vis2'],
self.poldict['vis3'],self.poldict['vis4']]:
if self.polrep=='stokes':
data[f] = np.hstack((self.data[f], np.conj(self.data[f])))
elif self.polrep=='circ':
if f in ['rrvis','llvis']:
data[f] = np.hstack((self.data[f], np.conj(self.data[f])))
elif f=='rlvis':
data[f] = np.hstack((self.data['rlvis'], np.conj(self.data['lrvis'])))
elif f=='lrvis':
data[f] = np.hstack((self.data['lrvis'], np.conj(self.data['rlvis'])))
else:
raise Exception("polrep must be either 'stokes' or 'circ'")
else:
data[f] = np.hstack((self.data[f], self.data[f]))
# Sort the data by time
data = data[np.argsort(data['time'])]
return data
def tlist(self, conj=False, t_gather=0., scan_gather=False):
"""Group the data in a list of equal time observation datatables.
Args:
conj (bool): True if tlist_out includes conjugate baselines.
t_gather (float): Grouping timescale (in seconds). 0.0 indicates no grouping.
scan_gather (bool): If true, gather data into scans
Returns:
(list): a list of data tables containing time-partitioned data
"""
if conj:
data = self.data_conj()
else:
data = self.data
# partition the data by time
datalist = []
if t_gather <= 0.0 and not scan_gather:
# Only group measurements at the same time
for key, group in it.groupby(data, lambda x: x['time']):
datalist.append(np.array([obs for obs in group]))
elif t_gather > 0.0 and not scan_gather:
# Group measurements in time
for key, group in it.groupby(data, lambda x: int(x['time']/(t_gather/3600.0))):
datalist.append(np.array([obs for obs in group]))
else:
# Group measurements by scan
if np.any(self.scans == None) or len(self.scans) == 0:
print("No scan table in observation. Adding scan table before gathering...")
self.add_scans()
for key, group in it.groupby(data, lambda x: np.searchsorted(self.scans[:,0],x['time'])):
datalist.append(np.array([obs for obs in group]))
return np.array(datalist)
def bllist(self, conj=False):
"""Group the data in a list of same baseline datatables.
Args:
conj (bool): True if tlist_out includes conjugate baselines.
Returns:
(list): a list of data tables containing baseline-partitioned data
"""
if conj:
data = self.data_conj()
else:
data = self.data
# partition the data by baseline
datalist = []
idx = np.lexsort((data['t2'], data['t1']))
for key, group in it.groupby(data[idx], lambda x: set((x['t1'], x['t2']))):
datalist.append(np.array([obs for obs in group]))
return np.array(datalist)
def unpack_bl(self, site1, site2, fields, ang_unit='deg', debias=False, timetype=False):
"""Unpack the data over time on the selected baseline site1-site2.
Args:
site1 (str): First site name
site2 (str): Second site name
fields (list): list of unpacked quantities from available quantities in FIELDS
ang_unit (str): 'deg' for degrees and 'rad' for radian phases
debias (bool): True to debias visibility amplitudes
timetype (str): 'GMST' or 'UTC'
Returns:
(numpy.recarray): unpacked numpy array with data in fields requested
"""
if timetype==False:
timetype=self.timetype
# If we only specify one field
if timetype not in ['GMST','UTC','utc','gmst']:
raise Exception("timetype should be 'GMST' or 'UTC'!")
allfields = ['time']
if not isinstance(fields, list): allfields.append(fields)
else:
for i in range(len(fields)): allfields.append(fields[i])
# Get the data from data table on the selected baseline
allout = []
tlist = self.tlist(conj=True)
for scan in tlist:
for obs in scan:
if (obs['t1'], obs['t2']) == (site1, site2):
obs = np.array([obs])
out = self.unpack_dat(obs, allfields, ang_unit=ang_unit, debias=debias, timetype=timetype)
allout.append(out)
#if len(allout)==0:
#raise Exception("Baseline %s-%s has no data!"%(site1,site2))
return np.array(allout)
def unpack(self, fields, mode='all', ang_unit='deg', debias=False, conj=False, timetype=False):
"""Unpack the data for the whole observation .
Args:
fields (list): list of unpacked quantities from availalbe quantities in FIELDS
mode (str): 'all' returns all data in single table, 'time' groups output by equal time, 'bl' groups by baseline
ang_unit (str): 'deg' for degrees and 'rad' for radian phases
debias (bool): True to debias visibility amplitudes
conj (bool): True to include conjugate baselines
timetype (str): 'GMST' or 'UTC'
Returns:
(numpy.recarray): unpacked numpy array with data in fields requested
"""
if not mode in ('time', 'all', 'bl'):
raise Exception("possible options for mode are 'time', 'all' and 'bl'")
# If we only specify one field
if not isinstance(fields, list): fields = [fields]
if mode=='all':
if conj:
data = self.data_conj()
else:
data = self.data
allout=self.unpack_dat(data, fields, ang_unit=ang_unit, debias=debias,timetype=timetype)
elif mode=='time':
allout=[]
tlist = self.tlist(conj=True)
for scan in tlist:
out=self.unpack_dat(scan, fields, ang_unit=ang_unit, debias=debias,timetype=timetype)
allout.append(out)
elif mode=='bl':
allout = []
bllist = self.bllist()
for bl in bllist:
out = self.unpack_dat(bl, fields, ang_unit=ang_unit, debias=debias,timetype=timetype)
allout.append(out)
return allout
def unpack_dat(self, data, fields, conj=False, ang_unit='deg', debias=False, timetype=False):
"""Unpack the data in a passed data recarray.
Args:
data (numpy.recarray): data recarray of format DTPOL_STOKES or DTPOL_CIRC
fields (list): list of unpacked quantities from availalbe quantities in FIELDS
conj (bool): True to include conjugate baselines
ang_unit (str): 'deg' for degrees and 'rad' for radian phases
debias (bool): True to debias visibility amplitudes
timetype (str): 'GMST' or 'UTC'
Returns:
(numpy.recarray): unpacked numpy array with data in fields requested
"""
if ang_unit=='deg': angle=DEGREE
else: angle = 1.0
# If we only specify one field
if type(fields) is str:
fields = [fields]
if not timetype:
timetype=self.timetype
if timetype not in ['GMST','UTC','gmst','utc']:
raise Exception("timetype should be 'GMST' or 'UTC'!")
# Get field data
allout = []
for field in fields:
if field in ["time","time_utc","time_gmst"]:
out = data['time']
ty='f8'
elif field in ["u","v","tint","tau1","tau2"]:
out = data[field]
ty = 'f8'
elif field in ["uvdist"]:
out = np.abs(data['u'] + 1j * data['v'])
ty = 'f8'
elif field in ["t1","el1","par_ang1","hr_ang1"]:
sites = data["t1"]
keys = [self.tkey[site] for site in sites]
tdata = self.tarr[keys]
out = sites
ty = 'U32'
elif field in ["t2","el2","par_ang2","hr_ang2"]:
sites = data["t2"]
keys = [self.tkey[site] for site in sites]
tdata = self.tarr[keys]
out = sites
ty = 'U32'
elif field in ['vis','amp','phase','snr','sigma','sigma_phase']:
ty = 'c16'
if self.polrep=='stokes':
out = data['vis']
sig = data['sigma']
elif self.polrep=='circ':
out = 0.5*(data['rrvis'] + data['llvis'])
sig = 0.5*np.sqrt(data['rrsigma']**2 + data['llsigma']**2)
elif field in ['qvis','qamp','qphase','qsnr','qsigma','qsigma_phase']:
ty = 'c16'
if self.polrep=='stokes':
out = data['qvis']
sig = data['qsigma']
elif self.polrep=='circ':
out = 0.5*(data['lrvis'] + data['rlvis'])
sig = 0.5*np.sqrt(data['lrsigma']**2 + data['rlsigma']**2)
elif field in ['uvis','uamp','uphase','usnr','usigma','usigma_phase']:
ty = 'c16'
if self.polrep=='stokes':
out = data['uvis']
sig = data['usigma']
elif self.polrep=='circ':
out = 0.5j*(data['lrvis'] - data['rlvis'])
sig = 0.5*np.sqrt(data['lrsigma']**2 + data['rlsigma']**2)
elif field in ['vvis','vamp','vphase','vsnr','vsigma','vsigma_phase']:
ty = 'c16'
if self.polrep=='stokes':
out = data['vvis']
sig = data['vsigma']
elif self.polrep=='circ':
out = 0.5*(data['rrvis'] - data['llvis'])
sig = 0.5*np.sqrt(data['rrsigma']**2 + data['llsigma']**2)
elif field in ['pvis','pamp','pphase','psnr','psigma','psigma_phase']:
ty = 'c16'
if self.polrep=='stokes':
out = data['qvis'] + 1j * data['uvis']
sig = np.sqrt(data['qsigma']**2 + data['usigma']**2)
elif self.polrep=='circ':
out = data['rlvis']
sig = data['rlsigma']
elif field in ['m','mamp','mphase','msnr','msigma','msigma_phase']:
ty = 'c16'
if self.polrep=='stokes':
out = (data['qvis'] + 1j * data['uvis'])/data['vis']
sig = merr(data['sigma'], data['qsigma'], data['usigma'], data['vis'], out)
elif self.polrep=='circ':
out = 2 * data['rlvis'] / (data['rrvis'] + data['llvis'])
sig = merr2(data['rlsigma'], data['rrsigma'], data['llsigma'], 0.5*(data['rrvis']+data['llvis']), out)
elif field in ['rrvis', 'rramp', 'rrphase', 'rrsnr', 'rrsigma', 'rrsigma_phase']:
ty = 'c16'
if self.polrep=='stokes':
out = data['vis'] + data['vvis']
sig = np.sqrt(data['sigma']**2 + data['vsigma']**2)
elif self.polrep=='circ':
out = data['rrvis']
sig = data['rrsigma']
elif field in ['llvis', 'llamp', 'llphase', 'llsnr', 'llsigma', 'llsigma_phase']:
ty = 'c16'
if self.polrep=='stokes':
out = data['vis'] - data['vvis']
sig = np.sqrt(data['sigma']**2 + data['vsigma']**2)
elif self.polrep=='circ':
out = data['llvis']
sig = data['llsigma']
elif field in ['rlvis', 'rlamp', 'rlphase', 'rlsnr', 'rlsigma', 'rlsigma_phase']:
ty = 'c16'
if self.polrep=='stokes':
out = data['qvis'] + 1j*data['uvis']
sig = np.sqrt(data['qsigma']**2 + data['usigma']**2)
elif self.polrep=='circ':
out = data['rlvis']
sig = data['rlsigma']
elif field in ['lrvis', 'lramp', 'lrphase', 'lrsnr', 'lrsigma', 'lrsigma_phase']:
ty = 'c16'
if self.polrep=='stokes':
out = data['qvis'] - 1j*data['uvis']
sig = np.sqrt(data['qsigma']**2 + data['usigma']**2)
elif self.polrep=='circ':
out = data['lrvis']
sig = data['lrsigma']
else: raise Exception("%s is not a valid field \n" % field +
"valid field values are: " + ' '.join(FIELDS))
if field in ["time_utc"] and timetype=='GMST':
out = gmst_to_utc(out, self.mjd)
if field in ["time_gmst"] and timetype=='UTC':
out = utc_to_gmst(out, self.mjd)
# Compute elevation and parallactic angles
if field in ["el1","el2","hr_ang1","hr_ang2","par_ang1","par_ang2"]:
if self.timetype=='GMST':
times_sid = data['time']
else:
times_sid = utc_to_gmst(data['time'], self.mjd)
thetas = np.mod((times_sid - self.ra)*HOUR, 2*np.pi)
coords = recarr_to_ndarr(tdata[['x','y','z']],'f8')
el_angle = elev(earthrot(coords, thetas), self.sourcevec())
latlon = xyz_2_latlong(coords)
hr_angles = hr_angle(times_sid*HOUR, latlon[:,1], self.ra*HOUR)
if field in ["el1","el2"]:
out = el_angle/angle
ty = 'f8'
if field in ["hr_ang1","hr_ang2"]:
out = hr_angles/angle
ty = 'f8'
if field in ["par_ang1","par_ang2"]:
par_ang = par_angle(hr_angles, latlon[:,0], self.dec*DEGREE)
out = par_ang/angle
ty = 'f8'
# Get arg/amps/snr
if field in ["amp", "qamp", "uamp","vamp","pamp","mamp","rramp","llamp","rlamp","lramp"]:
out = np.abs(out)
if debias:
out = amp_debias(out, sig)
ty = 'f8'
elif field in ["sigma","qsigma","usigma","vsigma","psigma","msigma",
"rrsigma","llsigma","rlsigma","lrsigma"]:
out = np.abs(sig)
ty = 'f8'
elif field in ["phase", "qphase", "uphase", "vphase","pphase",
"mphase","rrphase","llphase","lrphase","rlphase"]:
out = np.angle(out)/angle
ty = 'f8'
elif field in ["sigma_phase","qsigma_phase","usigma_phase",
"vsigma_phase","psigma_phase","msigma_phase",
"rrsigma_phase","llsigma_phase","rlsigma_phase","lrsigma_phase"]:
out = np.abs(sig)/np.abs(out)/angle
ty = 'f8'
elif field in ["snr", "qsnr", "usnr", "vsnr", "psnr", "msnr","rrsnr","llsnr","rlsnr","lrsnr"]:
out = np.abs(out)/np.abs(sig)
ty = 'f8'
# Reshape and stack with other fields
out = np.array(out, dtype=[(field, ty)])
if len(allout) > 0:
allout = rec.merge_arrays((allout, out), asrecarray=True, flatten=True)
else:
allout = out
return allout
def sourcevec(self):
"""Return the source position vector in geocentric coordinates at 0h GMST.
Args:
Returns:
(numpy.array): normal vector pointing to source in geocentric coordinates (m)
"""
sourcevec = np.array([np.cos(self.dec*DEGREE), 0, np.sin(self.dec*DEGREE)])
return sourcevec
def res(self):
"""Return the nominal resolution (1/longest baseline) of the observation in radians.
Args:
Returns:
(float): normal array resolution in radians
"""
res = 1.0/np.max(self.unpack('uvdist')['uvdist'])
return res
def split_obs(self):
"""Split single observation into multiple observation files, one per scan.
Args:
Returns:
(list): list of single-scan Obsdata objects
"""
print("Splitting Observation File into " + str(len(self.tlist())) + " scans")
arglist, argdict = self.obsdata_args()
# note that the tarr of the output includes all sites, even those that don't participate in the scan
splitlist = []
for tdata in self.tlist():
arglist[DATPOS] = tdata
splitlist.append(Obsdata(*arglist, **argdict))
return splitlist
def chisq(self, im, dtype='vis', pol='I', ttype='nfft', mask=[], **kwargs):
#cp_uv_min=False,
#debias=True, systematic_noise=0.0, systematic_cphase_noise=0.0, maxset=False,
#ttype='nfft',fft_pad_factor=2):
"""Give the reduced chi^2 of the observation for the specified image and datatype.
Args:
im (Image): image to test chi^2
dtype (str): data type of chi^2 (e.g., 'vis', 'amp', 'bs', 'cphase')
pol (str): polarization type ('I', 'Q', 'U', 'V', 'LL', 'RR', 'LR', or 'RL'
mask (arr): mask of same dimension as im.imvec to screen out pixels in chi^2 computation
ttype (str): if "fast" or "nfft" use FFT to produce visibilities. Else "direct" for DTFT
fft_pad_factor (float): zero pad the image to fft_pad_factor * image size in FFT
conv_func ('str'): The convolving function for gridding; options are 'gaussian', 'pill', and 'cubic'
p_rad (int): The pixel radius for the convolving function in gridding for FFTs
order ('str'): Interpolation order for sampling the FFT
systematic_noise (float): a fractional systematic noise tolerance to add to thermal sigmas
snrcut (float): a snr cutoff for including data in the chi^2 sum
debias (bool): if True then apply debiasing to amplitudes/closure amplitudes
weighting (str): 'natural' or 'uniform'
systematic_cphase_noise (float): a value in degrees to add to the closure phase sigmas
cp_uv_min (float): flag baselines shorter than this before forming closure quantities
maxset (bool): if True, use maximal set instead of minimal for closure quantities
Returns:
(float): image chi^2
"""
# TODO -- should import this at top, but the circular dependencies create a mess...
import ehtim.imaging.imager_utils as iu
if pol not in im._imdict.keys():
raise Exception(pol + ' is not in the current image. Consider changing the polarization basis of the image.')
(data, sigma, A) = iu.chisqdata(self, im, mask, dtype, pol=pol, ttype=ttype, **kwargs)
chisq = iu.chisq(im._imdict[pol], A, data, sigma, dtype, ttype=ttype, mask=mask)
return chisq
def recompute_uv(self):
"""Recompute u,v points using observation times and metadata
Args:
Returns:
(Obsdata): New Obsdata object containing the same data with recomputed u,v points
"""
times = self.data['time']
site1 = self.data['t1']
site2 = self.data['t2']
arr = ehtim.array.Array(self.tarr)
print ("Recomputing U,V Points using MJD %d \n RA %e \n DEC %e \n RF %e GHz"
% (self.mjd, self.ra, self.dec, self.rf/1.e9))
(timesout,uout,vout) = compute_uv_coordinates(arr, site1, site2, times,
self.mjd, self.ra, self.dec,self.rf,
timetype=self.timetype, elevmin=0, elevmax=90)
if len(timesout) != len(times):
raise Exception("len(timesout) != len(times) in recompute_uv: check elevation limits!!")
datatable = self.data.copy()
datatable['u'] = uout
datatable['v'] = vout
arglist, argdict = self.obsdata_args()
arglist[DATPOS] = np.array(datatable)
out = Obsdata(*arglist, **argdict)
return out
def avg_coherent(self, inttime, scan_avg=False):
"""Coherently average data along u,v tracks in chunks of length inttime (sec)
Args:
inttime (float): coherent integration time in seconds
scan_avg (bool): if True, average over scans in self.scans instead of intime
Returns:
(Obsdata): Obsdata object containing averaged data
"""
if (scan_avg==True)&(getattr(self.scans, "shape", None) is None or len(self.scans) == 0):
print('No scan data, ignoring scan_avg!')
scan_avg=False
if inttime <= 0.0 and scan_avg == False:
print('No averaging done!')
return self.copy()
vis_avg = coh_avg_vis(self,dt=inttime,return_type='rec',err_type='predicted',scan_avg=scan_avg)
arglist, argdict = self.obsdata_args()
arglist[DATPOS] = vis_avg
out = Obsdata(*arglist, **argdict)
return out
def avg_incoherent(self, inttime, scan_avg=False, debias=True, err_type='predicted'):
"""Incoherently average data along u,v tracks in chunks of length inttime (sec)
Args:
inttime (float): incoherent integration time in seconds
scan_avg (bool): if True, average over scans in self.scans instead of intime
debias (bool): if True, debias the averaged amplitudes
err_type (str): 'predicted' or 'measured'
Returns:
(Obsdata): Obsdata object containing averaged data
"""
print('Incoherently averaging data, putting phases to zero!')
amp_rec = incoh_avg_vis(self,dt=inttime,debias=debias,scan_avg=scan_avg,return_type='rec',rec_type='vis',err_type=err_type)
arglist, argdict = self.obsdata_args()
arglist[DATPOS] = amp_rec
out = Obsdata(*arglist, **argdict)
return out
def add_amp(self, avg_time=0, scan_avg=False, debias=True, err_type='predicted', return_type='rec', round_s=0.1, snrcut=0.):
"""Adds attribute self.amp: aan amplitude table with incoherently averaged amplitudes
Args:
avg_time (float): incoherent integration time in seconds
scan_avg (bool): if True, average over scans in self.scans instead of intime
debias (bool): if True then apply debiasing
err_type (str): 'predicted' or 'measured'
return_type: data frame ('df') or recarray ('rec')
round_s (float): accuracy of datetime object in seconds
snrcut (float): flag amplitudes with snr lower than this
"""
# Get spacing between datapoints in seconds
if len(set([x[0] for x in list(self.unpack('time'))])) > 1:
tint0 = np.min(np.diff(np.asarray(sorted(list(set([x[0] for x in list(self.unpack('time'))]))))))*3600.
else:
tint0 = 0.0
if avg_time <= tint0:
adf = make_amp(self, debias=debias, round_s=round_s)
if return_type=='rec':
adf = df_to_rec(adf,'amp')
print("Updated self.amp: no averaging")