-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathpolconvert_standalone.py
executable file
·2927 lines (2504 loc) · 114 KB
/
polconvert_standalone.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
# Copyright (c) Ivan Marti-Vidal 2012-2023
# EU ALMA Regional Center. Nordic node.
# Universitat de Valencia (Spain)
#
# 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 2 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/>,
# or write to the Free Software Foundation, Inc.,
# 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
#
#
# Redistribution and use in source and binary forms, with or without
# modification, are permitted provided that the following conditions are
# met:
#
# a. Redistributions of source code must retain the above copyright
# notice, this list of conditions and the following disclaimer.
# b. Redistributions in binary form must reproduce the above copyright
# notice, this list of conditions and the following disclaimer in the
# documentation and/or other materials provided with the
# distribution.
# c. Neither the name of the author nor the names of contributors may
# be used to endorse or promote products derived from this software
# without specific prior written permission.
#
#
# THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
# "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
# LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
# A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
# OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
# SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
# LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
# DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
# THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
# (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
# OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
#
#
#
from __future__ import absolute_import
from __future__ import print_function
__version__ = "2.0.7 " # 7 characters
date = "Aug 14, 2023"
################
# Import all necessary modules.
########
# Execute twice, to avoid the silly (and harmless)
# error regarding the different API versions of
# numpy between the local system and CASA:
import os, sys
try:
mypath = os.path.dirname(__file__)
sys.path.append(mypath)
import _PolConvert as PC
goodclib = True
print("\nC++ shared library loaded successfully (first try)\n")
except:
goodclib = False
print("\n There has been an error related to the numpy")
print(" API version used by CASA. This is related to PolConvert")
print(" (which uses the API version of your system) and should")
print(" be *harmless*.\n")
if not goodclib:
try:
import _PolConvert as PC
goodclib = True
print("\nC++ shared library loaded successfully (2nd try)\n")
except:
goodclib = False
print("\nNo, the shared library did not load successfully\n")
#############
if not goodclib:
try:
import _PolConvert as PC
goodclib = True
print("\nC++ shared library loaded successfully (3rd try)\n")
except:
goodclib = False
print("\nNo, the shared library still did not load successfully\n")
#############
### ... and we could keep on trying! (with an infinite loop! XD ).
import os, sys, shutil, re
import gc
import time
import struct as stk
import scipy.optimize as spopt
import scipy.interpolate as spint
import numpy as np
import pylab as pl
import datetime as dt
import sys
import pickle as pk
def polconvert(
IDI="",
OUTPUTIDI="",
DiFXinput="",
DiFXcalc="",
doIF=[],
linAntIdx=[],
Range=[],
XYadd={},
XYdel={},
XYratio={},
usePcal={},
swapXY=[],
swapRL=False,
feedRotation=[],
correctParangle=False,
IDI_conjugated=False,
plotIF=[],
plotRange=[],
plotAnt="",
excludeAnts=[],
excludeBaselines=[],
doSolve=-1,
solint=[1, 1],
doTest=True,
npix=-1,
solveAmp=True,
solveMethod="COBYLA",
calstokes=[1.0, 0.0, 0.0, 0.0],
calfield=-1,
plotSuffix="",
pcalSuffix="",
ALMAstuff=[],
saveArgs=False,
amp_norm=0.0,
IFoffset=0,
AC_MedianWindow=0,
XYpcalMode="bandpass",
UVTaper=1.e9,
useRates = False,
useDelays = False,
mounts = {}
):
"""POLCONVERT - STANDALONE VERSION 2.0.7
Similar parameters as the method defined in polconvert_CASA. The parameters specific
of this task (i.e., polconvert_standalone::polconvert) may be useful to parallelize polconvert:
plotSuffix: Suffix to add to the names of the plot files.
pcalSuffix: Suffix to add to the names of the original pcal files in the difx directories.
ALMAstuff: List of ALMA-specific parameters passed by polconvert_CASA::polconvert.
Other new parameters are:
IFoffset: treats the IF with index "i" as the same IF with index "i+IFoffset". This is
useful to process correctly the autocorrelations with multi-file has been
used in the observations.
AC_MedianWindow: Size (in channels) of a median window filter to be applied to the
autocorrelations. This is useful if we want to filter out phasecal peaks
(for a better fringe normalization).
XYpcalModel: Controls how the X-Y pcal differences are interpolated in frequency.
It can be either "bandpass" (default) or "multitone" (similar to fourfit's).
UVTaper: Apply a UVtaper (in meters) to the X-pol gain estimator (to avoid too high
cross-polarization fringes).
useRates: If True, the pol-wise weighted rate average will be use to increase
(hopefully) the coherence time of the data when the x-pol gains are
being estimated.
useDelays: If True, the pol-wise delay average will be used to increase (hopefully)
the SNR over wider frequency averages.
mounts: Antenna mounts (for SWIN case), given as a dictionary. The keywords are
the antenna codes and the elements are the mount types. The default antenna
mount is alt-azimuth. Recognized mounts are:
AZ (alt-az), EQ (equatorial), OB (orbit), XY (X-Y), NR (Nasmyth right), NL (Nasmyth left)
Example: if antenna Yebes-40m (code YB) has Nasmyth left, and all the other antennas
have alt-az mounts, then: mounts = {'YB':'NL'}
"""
if saveArgs:
ARGS = {
"IDI": IDI,
"OUTPUTIDI": OUTPUTIDI,
"DiFXinput": DiFXinput,
"DiFXcalc": DiFXcalc,
"doIF": doIF,
"linAntIdx": linAntIdx,
"Range": Range,
"XYadd": XYadd,
"XYdel": XYdel,
"XYratio": XYratio,
"usePcal": usePcal,
"swapXY": swapXY,
"swapRL": swapRL,
"feedRotation": feedRotation,
"correctParangle": correctParangle,
"IDI_conjugated": IDI_conjugated,
"plotIF": plotIF,
"plotRange": plotRange,
"plotAnt": plotAnt,
"excludeAnts": excludeAnts,
"excludeBaselines": excludeBaselines,
"doSolve": doSolve,
"solint": solint,
"doTest": doTest,
"npix": npix,
"amp_norm": amp_norm,
"solveAmp": solveAmp,
"solveMethod": solveMethod,
"calstokes": calstokes,
"calfield": calfield,
"ALMAstuff": ALMAstuff,
"IFoffset": IFoffset,
"AC_MedianWindow": AC_MedianWindow,
"XYpcalMode": XYpcalMode,
"UVTaper": UVTaper,
"useRates":useRates,
"useDelays":useDelays,
"mounts":mounts
}
OFF = open("PolConvert_standalone.last", "wb")
pk.dump(ARGS, OFF)
OFF.close()
## Set to non-ALMA:
if len(ALMAstuff) == 0:
PC.setPCMode(0)
# amp_norm=0.0
logName = "PolConvert%s.log" % plotSuffix
## Codes for antenna mounts:
MntCodes = {'AZ':0, 'EQ':1, 'OB':2, 'XY':3, 'NR':4, 'NL':5}
############################################
# this turns into the verbosity argument of _PolConvert.so
print("Entered polconvert_standalone::polconvert()")
DEBUG = False
if "POLCONVERTDEBUG" in os.environ:
if os.environ["POLCONVERTDEBUG"] == "True":
DEBUG = True
else:
DEBUG = False
print("DEBUG setting is " + str(DEBUG))
print("__name__ is " + __name__)
# Auxiliary function: print error and raise exception:
def printError(msg):
print(msg, "\n")
lfile = open(logName, "a")
print("\n" + msg + "\n", file=lfile)
lfile.close()
sys.stdout.flush()
raise Exception(msg)
# Auxiliary function: print message (terminal and log):
def printMsg(msg, doterm=True, dolog=True):
if doterm:
print(msg)
if dolog:
lfile = open(logName, "a")
print(msg, file=lfile)
lfile.close()
# Auxiliary function: derive job label from DiFXinput
def jobLabel(inputname):
label = inputname
try:
label = re.sub(".input", "", os.path.basename(label))
except:
pass
return label
tic = time.time()
greet = """
##########################################################################'
# POLCONVERT -- version. #'
# Please, add the POLCONVERT reference to your publications: #'
# #'
# Marti-Vidal, Roy, Conway & Zensus 2016, A&A, 587, 143 #'
# #'
##########################################################################'
"""
greetings = re.sub("version", __version__, greet)
printMsg(greetings, dolog=False)
printMsg("\n\nPOLCONVERT - VERSION %s" % __version__, doterm=False)
#########################################
# DO SOME SANITY CHECKS OF PARAMETERS
try:
doSolve = float(doSolve)
except:
printError("ERROR! doSolve should be a float!")
scipyMethods = ["COBYLA", "Nelder-Mead", "Powell", "BFGS", "Newton-CG", "SLSQP"]
allMethods = scipyMethods + ["gradient", "Levenberg-Marquardt"]
if solveMethod not in allMethods:
printError("ERROR! 'solveMethod' must be any of: %s" % (", ".join(allMethods)))
if type(calstokes) is not list:
printError("ERROR! Wrong calstokes!")
elif len(calstokes) != 4 and doSolve > 0.0:
printError(
"ERROR! calstokes should have 4 elements, not %d (%s)!"
% (len(calstokes), str(calstokes))
)
for ii, item in enumerate(calstokes):
try:
calstokes[ii] = float(item)
except:
printError(
"ERROR! calstokes should only have float elements; got %s!"
% str(type(item))
)
Stokes = list(calstokes)
if doSolve >= 0.0 and (
Stokes[0] <= 0.0
or Stokes[0] < np.sqrt(Stokes[1] ** 2.0 + Stokes[2] ** 2.0 + Stokes[3] ** 2.0)
):
printError("ERROR! Inconsistent Stokes parameters!")
# Will implement solveQU soon!
solveQU = False
# calfield = -1
if type(solveQU) is not bool:
printError("ERROR! Wrong solveQU!")
if type(calfield) is not int:
printError("ERROR! Wrong calfield!")
doConj = True
if type(IDI_conjugated) is not bool:
printError("ERROR! IDI_cojugated should be a boolean!")
else:
doConj = IDI_conjugated
if type(swapRL) is not bool:
printError("ERROR! swapRL should be a boolean!")
if type(doIF) is not list:
printError("ERROR! doIF should be a list of integers!")
else:
for elem in doIF:
if type(elem) is not int:
printError("ERROR! doIF should be a list of integers!")
if elem == 0:
printError(
"ERROR! IF numbers are given in AIPS (and FITS-IDI) convention\n(i.e., starting from 1; not from 0).\n"
)
if type(doTest) is not bool:
printError("ERROR! doTest should be a boolean!")
if doTest:
printMsg("Will only compute, but not update the output file(s)", doterm=False)
if type(linAntIdx) is not list:
printError("ERROR! linAntIdx should be a list of antenna indices or names!")
for elem in linAntIdx:
if type(elem) not in [int, str]:
printError("ERROR! linAntIdx should be a list of antenna indices or names!")
if type(solint) is not list or (len(solint) not in [2, 3]):
printError(
"ERROR! solint (%s) must be a list of two/three numbers!" % str(solint)
)
else:
try:
solint[0] = int(solint[0])
solint[1] = int(solint[1])
except:
printError(
"ERROR! solint (%s) must be a list of two/three numbers!" % str(solint)
)
if len(solint) == 3:
solint[2] = float(solint[2])
else: # Default dt
solint.append(100.0)
nALMA = len(linAntIdx)
if not os.path.exists(IDI):
printError("ERROR! IDI file (or SWIN folder) does not exist!")
if type(OUTPUTIDI) is not str:
printError("ERROR! OUTPUTIDI should be a string!")
if type(plotIF) is int:
if plotIF > 0:
plotIF = [plotIF]
else:
plotIF = []
for pli in plotIF:
if type(pli) is not int:
printError("ERROR! plotIF should be an integer or a list of integers!")
for pli in plotIF:
if pli not in doIF:
printError("ERROR! Only converted IFs can be plotted!")
if type(usePcal) is not dict:
printError("Invalid format for usePcal! Should be a dictionary!\n")
pcant = 0
for k in usePcal.keys():
if type(usePcal[k]) is not bool:
printError("The elements of usePcal must be booleans!\n")
elif usePcal[k]:
pcant += 1
isPcalUsed = pcant > 0
if len(swapXY) == 0:
swapXY = [False for swp in range(nALMA)]
if len(swapXY) != nALMA:
printError(
"Invalid format for swapXY!\n Should be a list of booleans, as large as the number of linear-polarization VLBI stations!"
)
for sxy in swapXY:
if type(sxy) is not bool:
printError(
"Invalid format for swapXY!\n Should be a list of booleans, as large as the number of linear-polarization VLBI stations!"
)
if len(Range) not in [0, 8]:
printError(
"Invalid format for Range! Should be either an empty list or a list of 8 integers!"
)
for rr in Range:
if type(rr) is not int:
printError(
"Invalid format for Range! Should be either an empty list or a list of 8 integers!"
)
if len(plotRange) not in [0, 8]:
printError(
"Invalid format for Range! Should be either an empty list or a list of 8 integers!"
)
for rr in plotRange:
if type(rr) is not int:
printError(
"Invalid format for Range! Should be either an empty list or a list of 8 integers!"
)
#########################################
#######
# CHECK IF THIS IS A FITS-IDI FILE OR A SWIN DIR.:
if os.path.isdir(IDI):
isSWIN = True
printMsg("\n\nYou have asked to convert a set of SWIN files.")
if len(DiFXinput) == 0:
DiFXinput = IDI[:-4] + "input"
if len(DiFXcalc) == 0:
DiFXcalc = IDI[:-4] + "calc"
if not os.path.exists(DiFXinput) or not os.path.isfile(DiFXinput):
printError("Invalid DiFX input file! %s" % DiFXinput)
printMsg("Opening calc file... %s" % DiFXcalc)
try:
printMsg('Opening "%s"' % (DiFXcalc))
antcoords = []
calcsoucoords = [[], []]
antmounts = []
antcodes = []
calc = open(DiFXcalc)
lines = calc.readlines()
calc.close()
printMsg("Read %d lines from %s" % (len(lines), DiFXcalc))
for ii, line in enumerate(lines):
if "TELESCOPE" in line and "NAME" in line:
antcodes.append(line.split()[-1])
if "TELESCOPE" in line and "X (m):" in line:
antcoords.append(list(map(float, [ll.split()[-1] for ll in lines[ii : ii + 3]])))
printMsg(
"TELESCOPE %s AT X: %.3f ; Y: %.3f ; Z: %.3f"
% tuple([antcodes[-1]] + antcoords[-1])
)
# DEFAULT MOUNTS ARE ALT-AZ:
antmounts.append(0)
if "SOURCE" in line and " NAME: " in line:
SNAM = line.split()[-1]
calcsoucoords[0].append(float(lines[ii + 1].split()[-1]))
calcsoucoords[1].append(float(lines[ii + 2].split()[-1]))
printMsg(
"SOURCE %s AT RA: %.8f rad, DEC: %.8f rad"
% (SNAM, calcsoucoords[0][-1], calcsoucoords[1][-1])
)
antcoords = np.array(antcoords, dtype=np.float, order="C")
antmounts = np.array(antmounts, order="C")
calcsoucoords[0] = np.array(calcsoucoords[0], dtype=np.float, order="C")
calcsoucoords[1] = np.array(calcsoucoords[1], dtype=np.float, order="C")
printMsg("done parsing calc")
except Exception as ex:
printMsg(str(ex))
printMsg(
(
"WARNING! Invalid DiFX calc file '%s'!\n"
+ "PolConvert may not calibrate properly."
)
% DiFXcalc
)
if len(antmounts) == 0:
printError("ERROR! NO ANTENNAS FOUND IN CALC FILE!")
else:
printMsg("There are %i antennas." % len(antmounts))
## Read antenna mounts from dictionary:
for antName in mounts.keys():
if antName in antcodes:
if mounts[antName] not in MntCodes.keys():
printError("ERROR! Mount %s not recognized!"%mounts[antName])
else:
code = mounts[antName]
mnt = MntCodes[code]
antmounts[antcodes.index(antName)] = int(mnt)
printMsg("Setting mount %s (%i) to antenna %s"%(code,mnt,antName))
elif os.path.isfile(IDI):
isSWIN = False
# Currently, PCALs are not supported for FITS-IDI:
isPcalUsed = False
printMsg("\n\nYou have asked to convert a FITS-IDI file.")
printMsg("Reading array geometry...")
try:
from astropy.io import fits as pf
ffile = pf.open(IDI)
for ii, group in enumerate(ffile):
if group.name == "ARRAY_GEOMETRY":
grarr = ii
elif group.name == "SOURCE":
grsou = ii
raappUnit = (
ffile[grsou]
.data.columns[
[coln.name for coln in ffile[grsou].data.columns].index("RAAPP")
]
.unit
)
decappUnit = (
ffile[grsou]
.data.columns[
[coln.name for coln in ffile[grsou].data.columns].index("DECAPP")
]
.unit
)
soucoords = [
np.array(ffile[grsou].data["RAAPP"], dtype=np.float, order="C"),
np.array(ffile[grsou].data["DECAPP"], dtype=np.float, order="C"),
]
if raappUnit == "DEGREES":
soucoords[0] *= np.pi / 180.0
if decappUnit == "DEGREES":
soucoords[1] *= np.pi / 180.0
antcodes = [ff[:2] for ff in ffile["ANTENNA"].data["ANNAME"]]
print("ANTENNA NAMES: ")
print(",".join(antcodes))
ffile.close()
# THESE LINES FAIL IF ORBPARM IS PRESENT IN ARRAY GEOMETRY!
import _getAntInfo as gA
success = gA.getAntInfo(IDI)
if success != 0:
printError("ERROR GETTING FITS-IDI METADATA! ERR: %i" % success)
else:
antcoords = gA.getCoords()
antmounts = gA.getMounts()
# FIXME: need a way to find antenna codes if above fails:
# antcodes = ['%02i'%i for i in range(1,len(antmounts)+1)]
except:
printMsg(
"WARNING! This FITS-IDI file has missing information!\nPolConvert may not calibrate properly."
)
else:
printError("Invalid input data!")
######
######
# IF THIS IS A SWIN DATASET, READ THE INPUT FILE INTO
# A METADATA LIST:
if isSWIN:
printMsg("Reading the DiFX input file\n")
ifile = open(DiFXinput)
inputlines = ifile.readlines()
ifile.close()
FreqL = [inputlines.index(l) for l in inputlines if "FREQ TABLE" in l]
# ONLY ONE FREQ TABLE IS ALLOWED:
try:
fr = FreqL[0]
Nfreq = int(
list(filter(lambda x: "FREQ ENTRIES" in x, inputlines[fr + 1 :]))[
0
].split()[-1]
)
Nr = list(range(Nfreq))
except Exception as ex:
printMsg(str(ex))
printError("BAD input file!")
FrInfo = {
"FREQ (MHZ)": [0.0 for i in Nr],
"BW (MHZ)": [0.0 for i in Nr],
"SIDEBAND": ["U" for i in Nr],
"NUM CHANNELS": [1 for i in Nr],
"CHANS TO AVG": [1 for i in Nr],
"OVERSAMPLE FAC.": [1 for i in Nr],
"DECIMATION FAC.": [1 for i in Nr],
"SIGN": [1.0 for i in Nr],
}
# READ METADATA FOR ALL FREQUENCIES:
for entry in FrInfo.keys():
for line in inputlines[fr + 1 :]:
if entry in line:
index = int((line.split(":")[0]).split()[-1])
FrInfo[entry][index] = type(FrInfo[entry][0])(line.split(":")[-1])
# SORT OUT THE CHANNEL FREQUENCIES:
if len(doIF) == 0:
doIF = list(range(1, len(Nr) + 1))
metadata = []
IFchan = 0
for nu in Nr:
nu0 = FrInfo["FREQ (MHZ)"][nu]
bw = FrInfo["BW (MHZ)"][nu]
nchan = FrInfo["NUM CHANNELS"][nu]
# MAX. NUMBER OF CHANNELS:
chav = FrInfo["CHANS TO AVG"][nu]
if nu in doIF:
IFchan = max([IFchan, int(nchan / chav)])
sb = {True: 1.0, False: -1.0}[FrInfo["SIDEBAND"][nu] == "U"]
FrInfo["SIGN"][nu] = float(sb)
freqs = (
nu0 + np.linspace((sb-1.0)/2.0, (sb+1.0)/2.0, nchan//chav, endpoint=False)*bw)*1.0e6
if float(nchan // chav) != float(nchan / chav):
printMsg("linspace check chan: %d / %d = %f"%(nchan, chav, float(nchan / chav)))
freqs = (
nu0 + np.linspace((sb-1.0)/2.0, (sb+1.0)/2.0, nchan//chav, endpoint=False)*bw)*1.0e6
metadata.append(freqs)
#####
######
# IF THIS IS A FITS-IDI FILE, READ THE METADATA LIST:
else:
# READ FREQUENCY INFO:
from astropy.io import fits as pf
fitsf = pf.open(IDI)
nu0 = fitsf["FREQUENCY"].header["REF_FREQ"]
bw = fitsf["FREQUENCY"].header["CHAN_BW"]
nch = fitsf["FREQUENCY"].header[
"NO_CHAN"
] # *fitsf['FREQUENCY'].header['NO_BAND']
IFchan = nch
Nr = fitsf["FREQUENCY"].header["NO_BAND"]
sgn = {True: 1.0, False: -1.0}[bw > 0.0]
FrInfo = {"FREQ (MHZ)": [], "BW (MHZ)": [], "SIGN": [], "NUM CHANNELS": []}
if sgn:
FrInfo["SIDEBAND"] = ["U" for i in range(Nr)]
else:
FrInfo["SIDEBAND"] = ["L" for i in range(Nr)]
metadata = []
for i in range(Nr):
FrInfo["FREQ (MHZ)"] += [(nu0 + i * bw * nch) / 1.0e6]
FrInfo["BW (MHZ)"] += [bw * nch / 1.0e6]
FrInfo["SIGN"] += [sgn]
FrInfo["NUM CHANNELS"] += [int(nch)]
freqs = (
nu0 + np.linspace((sgn - 1.0) / 2.0, (sgn + 1.0) / 2.0, nch, endpoint=False)*bw)
metadata.append(freqs)
FrInfo["CHANS TO AVG"] = [1 for i in range(Nr)]
FrInfo["OVERSAMPLE FAC."] = [1 for i in range(Nr)]
FrInfo["DECIMATION FAC."] = [1 for i in range(Nr)]
if len(doIF) == 0:
doIF = list(range(1, 1 + fitsf["FREQUENCY"].header["NO_BAND"]))
fitsf.close()
# ANTENNAS TO PARTICIPATE IN THE GAIN ESTIMATES:
nTotAnt = len(antcoords)
calAnts = []
for exA in antcodes:
if exA not in excludeAnts:
calAnts.append(antcodes.index(exA) + 1)
else:
printMsg("Excluding antenna %s from the solution." % str(exA))
##z following section
if type(plotAnt) is str and len(plotAnt) == 0:
plotAnt = 1
try:
plotAnt = int(plotAnt)
except:
if plotAnt not in antcodes:
printError("Reference antenna %s is not found in metadata among %s!" % (str(plotAnt), antcodes))
else:
plotAnt = antcodes.index(plotAnt) + 1
if plotAnt in linAntIdx or antcodes[plotAnt - 1] in linAntIdx:
printMsg(
"WARNING: Plotting will involve autocorrelations. \nThis has not been fully tested!"
)
FlagBas1 = []
FlagBas2 = []
for fbi in excludeBaselines:
printMsg("Excluding baseline %s from solution." % str(fbi))
if fbi[0] in antcodes and fbi[1] in antcodes:
FlagBas1.append(
antcodes.index(fbi[0]) + 1
) ### = np.array([int(i[0]+1) for i in excludeBaselines])
FlagBas2.append(
antcodes.index(fbi[1]) + 1
) ### = np.array([int(i[1]+1) for i in excludeBaselines])
else:
printError(
"Unknown antenna(s) %s and/or %s in excludeBaselines!\n"
% (fbi[0], fbi[1])
)
FlagBas1 = np.array(FlagBas1, order="C", dtype=np.int32)
FlagBas2 = np.array(FlagBas2, order="C", dtype=np.int32)
if plotAnt not in calAnts:
if doSolve >= 0:
printError(
"ERROR! plotAnt/Reference antenna is NOT in list of calibratable antennas!"
)
else:
printMsg(
"plotAnt (%d) is not in antenna list, so plots will be missing"
% plotAnt
)
if type(feedRotation) is not list:
printError("feedRotation must be a list of numbers")
elif len(feedRotation) == 0:
feedRot = np.zeros(nTotAnt, order="C", dtype=np.float)
elif len(feedRotation) != nTotAnt:
printError("feedRotation must have %i entries!" % nTotAnt)
else:
feedRot = np.pi / 180.0 * np.array(feedRotation, dtype=np.float, order="C")
# Get the REAL number (and names) of linear-pol antennas in this dataset:
nALMATrue = 0
linAntIdxTrue = []
linAntNamTrue = []
OrigLinIdx = []
for i, idd in enumerate(linAntIdx):
if type(idd) is int and idd <= len(antcodes):
nALMATrue += 1
linAntIdxTrue.append(idd)
linAntNamTrue.append(antcodes[idd - 1])
OrigLinIdx.append(i)
elif idd in antcodes:
nALMATrue += 1
linAntNamTrue.append(idd)
linAntIdxTrue.append(antcodes.index(idd) + 1)
OrigLinIdx.append(i)
# Setting isLinear list of booleans:
isLinear = []
for antc in antcodes:
isLinear.append(antc in linAntNamTrue)
# Check that pcals are (or not) set for this antenna:
for idd in linAntNamTrue:
if idd not in usePcal.keys():
printMsg(
"WARNING: Antenna %s not in usePcal dictionary. Will NOT use Pcals for this antenna"
% idd
)
usePcal[idd] = False
printMsg("There are %i linear-polarization antennas in THIS dataset" % nALMATrue)
# COMPUTE TIME RANGES:
if len(plotRange) == 0:
plRan = np.array([0.0, 0.0])
plotFringe = False
else:
try:
plRan = np.array(
[
plotRange[0]
+ plotRange[1]/24.0
+ plotRange[2]/1440.0
+ plotRange[3]/86400.0,
plotRange[4]
+ plotRange[5]/24.0
+ plotRange[6]/1440.0
+ plotRange[7]/86400.0,
]
)
plotFringe = True
if len(plotIF) == 0:
plotIF = list(doIF)
except:
printError("Bad time range format for plotRange!")
if len(Range) == 0:
Ran = np.array([0.0, 1.0e20])
else:
try:
Ran = np.array(
[
Range[0] + Range[1] / 24.0 + Range[2] / 1440.0 + Range[3] / 86400.0,
Range[4] + Range[5] / 24.0 + Range[6] / 1440.0 + Range[7] / 86400.0,
]
)
except:
printError("Bad time range format for Range!")
#######
# WARNING! UNCOMMENT THIS IF NOT DEBUGGING!
if os.path.exists(OUTPUTIDI) and IDI != OUTPUTIDI:
printMsg("Will REMOVE the existing OUTPUT file (or directory)!\n")
printMsg("Copying IDI to OUTPUTIDI!\n")
os.system("rm -rf %s" % OUTPUTIDI)
os.system("cp -r %s %s" % (IDI, OUTPUTIDI))
elif not os.path.exists(OUTPUTIDI):
printMsg("Copying IDI to OUTPUTIDI!\n")
os.system("cp -r %s %s" % (IDI, OUTPUTIDI))
#
#######
if isSWIN:
OUTPUT_UNSORT = []
PHASECALS_UNSORT = []
walk = [f for f in os.walk(OUTPUTIDI)]
for subd in walk:
ADD2OUTPUT = [
os.path.join(subd[0], fi)
for fi in filter(lambda x: x.startswith("DIFX_"), subd[2])
]
OUTPUT_UNSORT += ADD2OUTPUT
for difxfile in ADD2OUTPUT:
prefix = (os.path.basename(difxfile).split(".")[0])[5:]
PHASECALS_UNSORT.append(
[
os.path.join(subd[0], fi)
for fi in filter(
lambda x: x.startswith("PCAL_%s" % prefix), subd[2]
)
]
)
if len(OUTPUT_UNSORT) == 0:
printError("No *.difx files found in directory!")
OUTPUT = []
PHASECALS = []
TOOUT = np.argsort(
[
float(os.path.basename(fi).split(".")[0].split("_")[-1])
for fi in OUTPUT_UNSORT
]
)
for i in range(len(OUTPUT_UNSORT)):
OUTPUT.append(OUTPUT_UNSORT[TOOUT[i]])
PHASECALS.append(PHASECALS_UNSORT[TOOUT[i]])
# Derive the days of observation of each difx file:
mjd = np.zeros(len(OUTPUT))
mjs = np.zeros(len(OUTPUT))
for i, fi in enumerate(OUTPUT):
mjd[i], mjs[i] = list(
map(float, ((os.path.basename(fi)).split(".")[0]).split("_")[1:3])
)
mjd0 = np.min(mjd)
mjp = Ran + mjd0
metadata.append(mjd0)
# Filter out files outside the computing time window:
t0 = mjd + mjs / 86400.0
i0 = np.logical_and(t0 <= mjp[1], t0 >= mjp[0])
OUTPUT = [OUTPUT[i] for i in range(len(i0)) if i0[i]]
# Get source coordinate for each file (if possible):
#CALCS = ["%s.calc" % os.path.dirname(ci)[:-5] for ci in OUTPUT] # NB the [:-5] fails when user suffix used like '.pcdifx' or '.difx_PC' like by EUVGOS_PY3/POLCONVERTER.py
CALCS = ["%s.calc" % os.path.dirname(ci)[:os.path.dirname(ci).rfind('.')] for ci in OUTPUT]
##OBSOLETE! NOW, WE ARE CONCATENATING SWIN METADATA TO ENSURE CONSISTENCY!
# There will be one source per SWIN file:
# soucoords = [np.ones(len(OUTPUT))*calcsoucoords[0][0], np.ones(len(OUTPUT))*calcsoucoords[0][0]]
# for sui,Fcalc in enumerate(CALCS):
# # print(sui,Fcalc)
# if os.path.exists(Fcalc):
# FCin = open(Fcalc)
# lines = FCin.readlines()
# FCin.close()
# for ii,line in enumerate(lines):
# if 'SOURCE' in line and ' NAME: ' in line:
# SNAM = line.split()[-1]
# soucoords[0][sui] = float(lines[ii+1].split()[-1])
# soucoords[1][sui] = float(lines[ii+2].split()[-1])
# printMsg('FOUND SOURCE %s AT RA: %.8f rad, DEC: %.8f rad'%(SNAM,soucoords[0][sui],soucoords[1][sui]))
CALC2OPEN = {True: CALCS[0], False: DiFXcalc}[len(ALMAstuff)==0]
FCin = open(CALC2OPEN)
lines = FCin.readlines()
FCin.close()
for line in lines:
if line.startswith("NUM SOURCES: "):
NSOU = int(line.split()[-1].replace("\n", ""))
soucoords = [np.ones(NSOU), np.ones(NSOU)]
break
for ii, line in enumerate(lines):
if "SOURCE" in line and " NAME: " in line:
tempLine = line.split()
SID = int(tempLine[1])
soucoords[0][SID] = float(lines[ii + 1].split()[-1])
soucoords[1][SID] = float(lines[ii + 2].split()[-1])
else:
metadata = []
OUTPUT = [OUTPUTIDI]
##########################################
#########
# COMPUTE XY delays:
XYdelF = [[0.0, 0.0] for i in range(nALMATrue)]
if type(XYdel) is not dict: # or len(XYdel) != nALMA:
printError(
"Invalid format for XYdel!\n"
) # Should be a LIST of numbers, as large as the number of linear-polarization VLBI stations!")
for i, doant in enumerate(linAntNamTrue):
if doant in XYdel.keys():