-
Notifications
You must be signed in to change notification settings - Fork 17
/
Copy pathnoisedbase.py
3276 lines (3202 loc) · 182 KB
/
noisedbase.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
# -*- coding: utf-8 -*-
"""
A python module for ambient noise data analysis based on ASDF database
:Methods:
aftan analysis (use pyaftan or aftanf77)
preparing data for surface wave tomography (Barmin's method, Eikonal/Helmholtz tomography)
stacking/rotation for cross-correlation results from ANXCorr and Seed2Cor
:Dependencies:
pyasdf and its dependencies
ObsPy and its dependencies
pyproj
Basemap
pyfftw 0.10.3 (optional)
:Copyright:
Author: Lili Feng
Graduate Research Assistant
CIEI, Department of Physics, University of Colorado Boulder
email: lili.feng@colorado.edu
"""
import pyasdf
import numpy as np
import matplotlib.pyplot as plt
import matplotlib.dates as mdates
import obspy
import warnings
import copy
import os, shutil
from functools import partial
import multiprocessing
import pyaftan
from subprocess import call
from obspy.clients.fdsn.client import Client
from mpl_toolkits.basemap import Basemap, shiftgrid, cm
import obspy.signal.array_analysis
from obspy.imaging.cm import obspy_sequential
import glob
from numba import jit, float32, int32, boolean, float64
import pyfftw
import time
sta_info_default = {'xcorr': 1, 'isnet': 0}
xcorr_header_default = {'netcode1': '', 'stacode1': '', 'netcode2': '', 'stacode2': '', 'chan1': '', 'chan2': '',
'npts': 12345, 'b': 12345, 'e': 12345, 'delta': 12345, 'dist': 12345, 'az': 12345, 'baz': 12345, 'stackday': 0}
xcorr_sacheader_default = {'knetwk': '', 'kstnm': '', 'kcmpnm': '', 'stla': 12345, 'stlo': 12345,
'kuser0': '', 'kevnm': '', 'evla': 12345, 'evlo': 12345, 'evdp': 0., 'dist': 0., 'az': 12345, 'baz': 12345,
'delta': 12345, 'npts': 12345, 'user0': 0, 'b': 12345, 'e': 12345}
monthdict = {1: 'JAN', 2: 'FEB', 3: 'MAR', 4: 'APR', 5: 'MAY', 6: 'JUN', 7: 'JUL', 8: 'AUG', 9: 'SEP', 10: 'OCT', 11: 'NOV', 12: 'DEC'}
@jit(float32[:](float32[:, :], float32[:, :], int32))
def _CalcRecCor(arr1, arr2, lagN):
"""compute the amplitude weight for the xcorr, used for amplitude correction
optimized by numba
==============================================================================
::: input parameters :::
arr1, arr2 - the input arrays from ft_*SAC_rec,
indicating holes in the original records
lagN - one-sided npts for xcorr
::: output :::
cor_rec - amplitude weight for the xcorr
==============================================================================
"""
N1 = arr1.shape[0]
N2 = arr2.shape[0]
# array indicating number of data points for xcorr, used for amplitude correction
cor_rec = np.zeros(int(2*lagN + 1), dtype=float)
for i in range(lagN+1):
cor_rec[lagN+i] = 0.
cor_rec[lagN-i] = 0.
for irec1 in range(N1):
for irec2 in range(N2):
if arr1[irec1, 0] >= arr2[irec2, 1] - i:
continue
if arr1[irec1, 1] <= arr2[irec2, 0] - i:
break
recB = max(arr1[irec1, 0], arr2[irec2, 0] - i)
recE = min(arr1[irec1, 1], arr2[irec2, 1] - i)
cor_rec[lagN+i] += recE - recB
for irec2 in range(N2):
if arr1[irec1, 0] >= arr2[irec2, 1] + i:
continue
if arr1[irec1, 1] <= arr2[irec2, 0] + i:
break
recB = max(arr1[irec1, 0], arr2[irec2, 0] + i)
recE = min(arr1[irec1, 1], arr2[irec2, 1] + i)
cor_rec[lagN-i] += recE - recB
cor_rec[lagN] /= 2.
return cor_rec
def _amp_ph_to_xcorr(amp1, amp2, ph1, ph2, sps = 1., lagtime = 3000.):
"""Convert amplitude and phase arrays to xcorr
==============================================================================
::: input parameters :::
amp1, ph1 - amplitude and phase data arrays for station(component) 1
amp2, ph2 - amplitude and phase data arrays for station(component) 2
sps - target sampling rate
lagtime - lag time for xcorr
::: output :::
out_data - xcorr
==============================================================================
"""
N = amp1.size
Ns = int(2*N - 1)
# cross-spectrum, conj(sac1)*(sac2)
x_sp = np.zeros(Ns, dtype=complex)
temp1 = np.zeros(N, dtype=complex)
temp2 = np.zeros(N, dtype=complex)
temp1.real = amp1*np.cos(ph1)
temp1.imag = amp1*np.sin(ph1)
temp2.real = amp2*np.cos(ph2)
temp2.imag = amp2*np.sin(ph2)
x_sp[:N] = temp2 * np.conj(temp1)
# perform inverse FFT with pyFFTW, much faster than numpy_fft, scipy.fftpack
out = pyfftw.interfaces.numpy_fft.ifft(x_sp)
seis_out = 2.*(out.real)
lagN = int(np.floor(lagtime*sps +0.5))
if lagN > Ns:
raise ValueError('Lagtime npts overflow!')
out_data = np.zeros(2*lagN+1, dtype=float)
out_data[lagN] = seis_out[0]
out_data[:lagN] = (seis_out[1:lagN+1])[::-1]
out_data[(lagN+1):] = (seis_out[Ns-lagN:])[::-1]
return out_data
def _amp_ph_to_xcorr_fast(amp1, amp2, ph1, ph2, fftw_plan, sps = 1., lagtime = 3000.):
"""Convert amplitude and phase arrays to xcorr
This is the fast version of _amp_ph_to_xcorr, a precomputed fftw_plan needs
to be prepared for speeding up
==============================================================================
::: input parameters :::
amp1, ph1 - amplitude and phase data arrays for station(component) 1
amp2, ph2 - amplitude and phase data arrays for station(component) 2
fftw_plan - precomputed fftw_plan
sps - target sampling rate
lagtime - lag time for xcorr
::: output :::
out_data - xcorr
==============================================================================
"""
N = amp1.size
Ns = int(2*N - 1)
# cross-spectrum, conj(sac1)*(sac2)
x_sp = np.zeros(Ns, dtype=complex)
out = np.zeros(Ns, dtype=complex)
temp1 = np.zeros(N, dtype=complex)
temp2 = np.zeros(N, dtype=complex)
temp1.real = amp1*np.cos(ph1)
temp1.imag = amp1*np.sin(ph1)
temp2.real = amp2*np.cos(ph2)
temp2.imag = amp2*np.sin(ph2)
x_sp[:N] = temp2 * np.conj(temp1)
# perform inverse FFT with pyFFTW, much faster than numpy_fft, scipy.fftpack
# the precomputed fftw_plan is used
fftw_plan.update_arrays(x_sp, out)
fftw_plan.execute()
seis_out = 2.*(out.real)
lagN = int(np.floor(lagtime*sps +0.5))
if lagN > Ns:
raise ValueError('Lagtime npts overflow!')
out_data = np.zeros(2*lagN+1, dtype=float)
out_data[lagN] = seis_out[0]
out_data[:lagN] = (seis_out[1:lagN+1])[::-1]
out_data[(lagN+1):] = (seis_out[Ns-lagN:])[::-1]
return out_data/Ns
class xcorr_pair(object):
""" An object to for ambient noise cross-correlation computation
=================================================================================================================
::: parameters :::
stacode1, netcode1 - station/network code for station 1
stacode2, netcode2 - station/network code for station 2
monthdir - month directory (e.g. 2019.JAN)
daylst - list includes the days for xcorr
=================================================================================================================
"""
def __init__(self, stacode1, netcode1, stacode2, netcode2, monthdir, daylst):
self.stacode1 = stacode1
self.netcode1 = netcode1
self.stacode2 = stacode2
self.netcode2 = netcode2
self.monthdir = monthdir
self.daylst = daylst
return
def print_info(self):
"""print the informations of this pair
"""
staid1 = self.netcode1 + '.' + self.stacode1
staid2 = self.netcode2 + '.' + self.stacode2
print '--- '+ staid1+'_'+staid2+' : '+self.monthdir+' '+str(len(self.daylst))+' days'
def convert_amph_to_xcorr(self, datadir, chans=['LHZ', 'LHE', 'LHN'], ftlen = True,\
tlen = 84000., mintlen = 20000., sps = 1., lagtime = 3000., CorOutflag = 0, \
fprcs = False, fastfft=True, verbose=False):
"""
Convert amplitude and phase files to xcorr
=================================================================================================================
::: input parameters :::
datadir - directory including data and output
chans - channel list
ftlen - turn (on/off) cross-correlation-time-length correction for amplitude
tlen - time length of daily records (in sec)
mintlen - allowed minimum time length for cross-correlation (takes effect only when ftlen = True)
sps - target sampling rate
lagtime - cross-correlation signal half length in sec
CorOutflag - 0 = only output monthly xcorr data, 1 = only daily, 2 or others = output both
fprcs - turn on/off (1/0) precursor signal checking, NOT implemented yet
fastfft - speeding up the computation by using precomputed fftw_plan or not
=================================================================================================================
"""
if verbose:
self.print_info()
staid1 = self.netcode1 + '.' + self.stacode1
staid2 = self.netcode2 + '.' + self.stacode2
month_dir = datadir+'/'+self.monthdir
monthly_xcorr = []
chan_size = len(chans)
init_common_header = False
xcorr_common_sacheader = xcorr_sacheader_default.copy()
lagN = np.floor(lagtime*sps +0.5) # npts for one-sided lag
stacked_day = 0
#---------------------------------------
# construct fftw_plan for speeding up
#---------------------------------------
if fastfft:
temp_pfx = month_dir+'/'+self.monthdir+'.'+str(self.daylst[0])+\
'/ft_'+self.monthdir+'.'+str(self.daylst[0])+'.'+staid1+'.'+chans[0]+'.SAC'
amp_ref = obspy.read(temp_pfx+'.am')[0]
Nref = amp_ref.data.size
Ns = int(2*Nref - 1)
temp_x_sp = np.zeros(Ns, dtype=complex)
temp_out = np.zeros(Ns, dtype=complex)
fftw_plan = pyfftw.FFTW(input_array=temp_x_sp, output_array=temp_out, direction='FFTW_BACKWARD',\
flags=('FFTW_MEASURE', ))
else:
Nref = 0
#-----------------
# loop over days
#-----------------
for day in self.daylst:
# input streams
st_amp1 = obspy.Stream()
st_ph1 = obspy.Stream()
st_amp2 = obspy.Stream()
st_ph2 = obspy.Stream()
# daily output streams
daily_xcorr = []
daydir = month_dir+'/'+self.monthdir+'.'+str(day)
# read amp/ph files
for chan in chans:
pfx1 = daydir+'/ft_'+self.monthdir+'.'+str(day)+'.'+staid1+'.'+chan+'.SAC'
pfx2 = daydir+'/ft_'+self.monthdir+'.'+str(day)+'.'+staid2+'.'+chan+'.SAC'
st_amp1 += obspy.read(pfx1+'.am')
st_ph1 += obspy.read(pfx1+'.ph')
st_amp2 += obspy.read(pfx2+'.am')
st_ph2 += obspy.read(pfx2+'.ph')
#-----------------------------
# define commone sac header
#-----------------------------
if not init_common_header:
tr1 = st_amp1[0]
tr2 = st_amp2[0]
xcorr_common_sacheader['kuser0'] = self.netcode1
xcorr_common_sacheader['kevnm'] = self.stacode1
xcorr_common_sacheader['knetwk'] = self.netcode2
xcorr_common_sacheader['kstnm'] = self.stacode2
# # # xcorr_common_sacheader['kcmpnm'] = chan1+chan2
xcorr_common_sacheader['evla'] = tr1.stats.sac.stla
xcorr_common_sacheader['evlo'] = tr1.stats.sac.stlo
xcorr_common_sacheader['stla'] = tr2.stats.sac.stla
xcorr_common_sacheader['stlo'] = tr2.stats.sac.stlo
dist, az, baz = obspy.geodetics.gps2dist_azimuth(tr1.stats.sac.stla, tr1.stats.sac.stlo,\
tr2.stats.sac.stla, tr2.stats.sac.stlo) # distance is in m
xcorr_common_sacheader['dist'] = dist/1000.
xcorr_common_sacheader['az'] = az
xcorr_common_sacheader['baz'] = baz
xcorr_common_sacheader['delta'] = 1./sps
xcorr_common_sacheader['npts'] = int(2*lagN + 1)
xcorr_common_sacheader['b'] = -float(lagN/sps)
xcorr_common_sacheader['e'] = float(lagN/sps)
xcorr_common_sacheader['user0'] = 1
init_common_header = True
skip_this_day = False
# compute cross-correlation
for ich1 in range(chan_size):
if skip_this_day:
break
for ich2 in range(chan_size):
# get data arrays
amp1 = st_amp1[ich1].data
ph1 = st_ph1[ich1].data
amp2 = st_amp2[ich2].data
ph2 = st_ph2[ich2].data
# quality control
if (np.isnan(amp1)).any() or (np.isnan(amp2)).any() or \
(np.isnan(ph1)).any() or (np.isnan(ph2)).any():
skip_this_day = True
break
if np.any(amp1 > 1e20) or np.any(amp2 > 1e20):
skip_this_day = True
break
# get amplitude correction array
Namp = amp1.size
if ftlen:
# npts for the length of the preprocessed daily record
Nrec = int(tlen*sps)
frec1 = daydir+'/ft_'+self.monthdir+'.'+str(day)+'.'+staid1+'.'+chans[ich1]+'.SAC_rec'
frec2 = daydir+'/ft_'+self.monthdir+'.'+str(day)+'.'+staid2+'.'+chans[ich1]+'.SAC_rec'
if os.path.isfile(frec1):
arr1= np.loadtxt(frec1)
if arr1.size == 2:
arr1 = arr1.reshape(1, 2)
else:
arr1= (np.array([0, Nrec])).reshape(1, 2)
if os.path.isfile(frec2):
arr2= np.loadtxt(frec2)
if arr2.size == 2:
arr2 = arr2.reshape(1, 2)
else:
arr2= (np.array([0, Nrec])).reshape(1, 2)
cor_rec = _CalcRecCor(arr1, arr2, np.int32(lagN))
# skip the day if the length of available data is too small
if cor_rec[0] < mintlen*sps or cor_rec[-1] < mintlen*sps:
skip_this_day = True
break
# skip the day if any data point has a weight of zero
if np.any(cor_rec == 0.):
skip_this_day = True
break
# comvert amp & ph files to xcorr
if fastfft and Namp == Nref:
out_data = _amp_ph_to_xcorr_fast(amp1=amp1, ph1=ph1, amp2=amp2, ph2=ph2, sps=sps,\
lagtime=lagtime, fftw_plan=fftw_plan)
else:
out_data = _amp_ph_to_xcorr(amp1=amp1, ph1=ph1, amp2=amp2, ph2=ph2, sps=sps, lagtime=lagtime)
# amplitude correction
if ftlen:
out_data /= cor_rec
out_data *= float(2*Namp - 1)
# end of computing individual xcorr
daily_xcorr.append(out_data)
# end loop over channels
if not skip_this_day:
if verbose:
print 'xcorr finished '+ staid1+'_'+staid2+' : '+self.monthdir+'.'+str(day)
# output daily xcorr
if CorOutflag != 0:
out_daily_dir = month_dir+'/COR_D/'+staid1
if not os.path.isdir(out_daily_dir):
os.makedirs(out_daily_dir)
for ich1 in range(chan_size):
for ich2 in range(chan_size):
i = 3*ich1 + ich2
out_daily_fname = out_daily_dir+'/COR_'+staid1+'_'+chans[ich1]+\
'_'+staid2+'_'+chans[ich2]+'_'+str(day)+'.SAC'
daily_header = xcorr_common_sacheader.copy()
daily_header['kcmpnm'] = chans[ich1]+chans[ich2]
sacTr = obspy.io.sac.sactrace.SACTrace(data = daily_xcorr[i], **daily_header)
sacTr.write(out_daily_fname)
# append to monthly data
if CorOutflag != 1:
# intilize
if stacked_day == 0:
for ich1 in range(chan_size):
for ich2 in range(chan_size):
i = 3*ich1 + ich2
monthly_xcorr.append(daily_xcorr[i])
# stacking
else:
for ich1 in range(chan_size):
for ich2 in range(chan_size):
i = 3*ich1 + ich2
monthly_xcorr[i] += daily_xcorr[i]
stacked_day += 1
# end loop over days
if CorOutflag != 1 and stacked_day != 0:
out_monthly_dir = month_dir+'/COR/'+staid1
if not os.path.isdir(out_monthly_dir):
try:
os.makedirs(out_monthly_dir)
except OSError:
i = 0
while(i < 10):
sleep_time = np.random.random()/10.
time.sleep(sleep_time)
if not os.path.isdir(out_monthly_dir):
try:
os.makedirs(out_monthly_dir)
break
except OSError:
pass
i += 1
for ich1 in range(chan_size):
for ich2 in range(chan_size):
i = 3*ich1 + ich2
out_monthly_fname = out_monthly_dir+'/COR_'+staid1+'_'+chans[ich1]+\
'_'+staid2+'_'+chans[ich2]+'.SAC'
monthly_header = xcorr_common_sacheader.copy()
monthly_header['kcmpnm'] = chans[ich1]+chans[ich2]
monthly_header['user0'] = stacked_day
sacTr = obspy.io.sac.sactrace.SACTrace(data = monthly_xcorr[i], **monthly_header)
sacTr.write(out_monthly_fname)
return
def amph_to_xcorr_for_mp(in_xcorr_pair, datadir, chans=['LHZ', 'LHE', 'LHN'], ftlen = True,\
tlen = 84000., mintlen = 20000., sps = 1., lagtime = 3000., CorOutflag = 0, \
fprcs = False, fastfft=True):
in_xcorr_pair.convert_amph_to_xcorr(datadir=datadir, chans=chans, ftlen = ftlen,\
tlen = tlen, mintlen = mintlen, sps = sps, lagtime = lagtime, CorOutflag = CorOutflag,\
fprcs = fprcs, fastfft=fastfft)
# # # in_xcorr_pair.print_info()
return
class beamforming_stream(obspy.Stream):
""" An object to for ambient noise cross-correlation computation
=================================================================================================================
::: parameters :::
stacode1, netcode1 - station/network code for station 1
stacode2, netcode2 - station/network code for station 2
monthdir - month directory (e.g. 2019.JAN)
daylst - list includes the days for xcorr
=================================================================================================================
"""
def get_addtional_info(self, datadir, bfUTCdate):
self.datadir = datadir
self.bfUTCdate = bfUTCdate
return
class noiseASDF(pyasdf.ASDFDataSet):
""" An object to for ambient noise cross-correlation analysis based on ASDF database
=================================================================================================================
version history:
Dec 6th, 2016 - first version
=================================================================================================================
"""
def print_info(self):
"""
print information of the dataset.
"""
outstr = '================================================= Ambient Noise Cross-correlation Database =================================================\n'
outstr += self.__str__()+'\n'
outstr += '--------------------------------------------------------------------------------------------------------------------------------------------\n'
if 'NoiseXcorr' in self.auxiliary_data.list():
outstr += 'NoiseXcorr - Cross-correlation seismogram\n'
if 'StaInfo' in self.auxiliary_data.list():
outstr += 'StaInfo - Auxiliary station information\n'
if 'DISPbasic1' in self.auxiliary_data.list():
outstr += 'DISPbasic1 - Basic dispersion curve, no jump correction\n'
if 'DISPbasic2' in self.auxiliary_data.list():
outstr += 'DISPbasic2 - Basic dispersion curve, with jump correction\n'
if 'DISPpmf1' in self.auxiliary_data.list():
outstr += 'DISPpmf1 - PMF dispersion curve, no jump correction\n'
if 'DISPpmf2' in self.auxiliary_data.list():
outstr += 'DISPpmf2 - PMF dispersion curve, with jump correction\n'
if 'DISPbasic1interp' in self.auxiliary_data.list():
outstr += 'DISPbasic1interp - Interpolated DISPbasic1\n'
if 'DISPbasic2interp' in self.auxiliary_data.list():
outstr += 'DISPbasic2interp - Interpolated DISPbasic2\n'
if 'DISPpmf1interp' in self.auxiliary_data.list():
outstr += 'DISPpmf1interp - Interpolated DISPpmf1\n'
if 'DISPpmf2interp' in self.auxiliary_data.list():
outstr += 'DISPpmf2interp - Interpolated DISPpmf2\n'
if 'FieldDISPbasic1interp' in self.auxiliary_data.list():
outstr += 'FieldDISPbasic1interp - Field data of DISPbasic1\n'
if 'FieldDISPbasic2interp' in self.auxiliary_data.list():
outstr += 'FieldDISPbasic2interp - Field data of DISPbasic2\n'
if 'FieldDISPpmf1interp' in self.auxiliary_data.list():
outstr += 'FieldDISPpmf1interp - Field data of DISPpmf1\n'
if 'FieldDISPpmf2interp' in self.auxiliary_data.list():
outstr += 'FieldDISPpmf2interp - Field data of DISPpmf2\n'
outstr += '============================================================================================================================================\n'
print outstr
return
def write_stationxml(self, staxml, source='CIEI'):
"""write obspy inventory to StationXML data file
"""
inv = obspy.core.inventory.inventory.Inventory(networks=[], source=source)
for staid in self.waveforms.list():
inv += self.waveforms[staid].StationXML
inv.write(staxml, format='stationxml')
return
def write_stationtxt(self, stafile):
"""write obspy inventory to txt station list(format used in ANXcorr and Seed2Cor)
"""
try:
auxiliary_info = self.auxiliary_data.StaInfo
isStaInfo = True
except:
isStaInfo = False
with open(stafile, 'w') as f:
for staid in self.waveforms.list():
stainv = self.waveforms[staid].StationXML
netcode = stainv.networks[0].code
stacode = stainv.networks[0].stations[0].code
lon = stainv.networks[0].stations[0].longitude
lat = stainv.networks[0].stations[0].latitude
if isStaInfo:
staid_aux = netcode+'/'+stacode
xcorrflag = auxiliary_info[staid_aux].parameters['xcorr']
f.writelines('%s %3.4f %3.4f %d %s\n' %(stacode, lon, lat, xcorrflag, netcode) )
else:
f.writelines('%s %3.4f %3.4f %s\n' %(stacode, lon, lat, netcode) )
return
def read_stationtxt(self, stafile, source='CIEI', chans=['LHZ', 'LHE', 'LHN'], dnetcode=None):
"""read txt station list
"""
sta_info = sta_info_default.copy()
with open(stafile, 'r') as f:
Sta = []
site = obspy.core.inventory.util.Site(name='01')
creation_date = obspy.core.utcdatetime.UTCDateTime(0)
inv = obspy.core.inventory.inventory.Inventory(networks=[], source=source)
total_number_of_channels = len(chans)
for lines in f.readlines():
lines = lines.split()
stacode = lines[0]
lon = float(lines[1])
lat = float(lines[2])
netcode = dnetcode
xcorrflag = None
if len(lines)==5:
try:
xcorrflag = int(lines[3])
netcode = lines[4]
except ValueError:
xcorrflag = int(lines[4])
netcode = lines[3]
if len(lines)==4:
try:
xcorrflag = int(lines[3])
except ValueError:
netcode = lines[3]
if netcode is None:
netsta = stacode
else:
netsta = netcode+'.'+stacode
if Sta.__contains__(netsta):
index = Sta.index(netsta)
if abs(self[index].lon-lon) >0.01 and abs(self[index].lat-lat) >0.01:
raise ValueError('incompatible station location:' + netsta+' in station list!')
else:
print 'WARNING: repeated station:' +netsta+' in station list!'
continue
channels = []
if lon>180.:
lon -= 360.
for chan in chans:
channel = obspy.core.inventory.channel.Channel(code=chan, location_code='01', latitude=lat, longitude=lon,
elevation=0.0, depth=0.0)
channels.append(channel)
station = obspy.core.inventory.station.Station(code=stacode, latitude=lat, longitude=lon, elevation=0.0,
site=site, channels=channels, total_number_of_channels = total_number_of_channels, creation_date = creation_date)
network = obspy.core.inventory.network.Network(code=netcode, stations=[station])
networks = [network]
inv += obspy.core.inventory.inventory.Inventory(networks=networks, source=source)
if netcode is None:
staid_aux = stacode
else:
staid_aux = netcode+'/'+stacode
if xcorrflag != None:
sta_info['xcorr'] = xcorrflag
self.add_auxiliary_data(data=np.array([]), data_type='StaInfo', path=staid_aux, parameters=sta_info)
print('Writing obspy inventory to ASDF dataset')
self.add_stationxml(inv)
print('End writing obspy inventory to ASDF dataset')
return
def read_stationtxt_ind(self, stafile, source='CIEI', chans=['LHZ', 'LHE', 'LHN'], s_ind=1, lon_ind=2, lat_ind=3, n_ind=0):
"""read txt station list, column index can be changed
"""
sta_info = sta_info_default.copy()
with open(stafile, 'r') as f:
Sta = []
site = obspy.core.inventory.util.Site(name='')
creation_date = obspy.core.utcdatetime.UTCDateTime(0)
inv = obspy.core.inventory.inventory.Inventory(networks=[], source=source)
total_number_of_channels= len(chans)
for lines in f.readlines():
lines = lines.split()
stacode = lines[s_ind]
lon = float(lines[lon_ind])
lat = float(lines[lat_ind])
netcode = lines[n_ind]
netsta = netcode+'.'+stacode
if Sta.__contains__(netsta):
index = Sta.index(netsta)
if abs(self[index].lon-lon) >0.01 and abs(self[index].lat-lat) >0.01:
raise ValueError('incompatible station location:' + netsta+' in station list!')
else:
print 'WARNING: repeated station:' +netsta+' in station list!'
continue
channels = []
if lon>180.:
lon -=360.
for chan in chans:
channel = obspy.core.inventory.channel.Channel(code=chan, location_code='01', latitude=lat, longitude=lon,
elevation=0.0, depth=0.0)
channels.append(channel)
station = obspy.core.inventory.station.Station(code=stacode, latitude=lat, longitude=lon, elevation=0.0,
site=site, channels=channels, total_number_of_channels = total_number_of_channels, creation_date = creation_date)
network = obspy.core.inventory.network.Network(code=netcode, stations=[station])
networks = [network]
inv += obspy.core.inventory.inventory.Inventory(networks=networks, source=source)
staid_aux = netcode+'/'+stacode
self.add_auxiliary_data(data=np.array([]), data_type='StaInfo', path=staid_aux, parameters=sta_info)
print('Writing obspy inventory to ASDF dataset')
self.add_stationxml(inv)
print('End writing obspy inventory to ASDF dataset')
return
def copy_stations(self, inasdffname, startdate=None, enddate=None, location=None, channel=None, includerestricted=False,
minlatitude=None, maxlatitude=None, minlongitude=None, maxlongitude=None, latitude=None, longitude=None, minradius=None, maxradius=None):
"""copy and renew station inventory given an input ASDF file
the function will copy the network and station names while renew other informations given new limitations
=======================================================================================================
::: input parameters :::
inasdffname - input ASDF file name
startdate, enddata - start/end date for searching
network - Select one or more network codes.
Can be SEED network codes or data center defined codes.
Multiple codes are comma-separated (e.g. "IU,TA").
station - Select one or more SEED station codes.
Multiple codes are comma-separated (e.g. "ANMO,PFO").
location - Select one or more SEED location identifiers.
Multiple identifiers are comma-separated (e.g. "00,01").
As a special case “--“ (two dashes) will be translated to a string of two space
characters to match blank location IDs.
channel - Select one or more SEED channel codes.
Multiple codes are comma-separated (e.g. "BHZ,HHZ").
minlatitude - Limit to events with a latitude larger than the specified minimum.
maxlatitude - Limit to events with a latitude smaller than the specified maximum.
minlongitude - Limit to events with a longitude larger than the specified minimum.
maxlongitude - Limit to events with a longitude smaller than the specified maximum.
latitude - Specify the latitude to be used for a radius search.
longitude - Specify the longitude to the used for a radius search.
minradius - Limit to events within the specified minimum number of degrees from the
geographic point defined by the latitude and longitude parameters.
maxradius - Limit to events within the specified maximum number of degrees from the
geographic point defined by the latitude and longitude parameters.
=======================================================================================================
"""
try:
starttime = obspy.core.utcdatetime.UTCDateTime(startdate)
except:
starttime = None
try:
endtime = obspy.core.utcdatetime.UTCDateTime(enddate)
except:
endtime = None
client = Client('IRIS')
init_flag = False
indset = pyasdf.ASDFDataSet(inasdffname)
for staid in indset.waveforms.list():
network = staid.split('.')[0]
station = staid.split('.')[1]
print 'Copying/renewing station inventory: '+ staid
if init_flag:
inv += client.get_stations(network=network, station=station, starttime=starttime, endtime=endtime, channel=channel,
minlatitude=minlatitude, maxlatitude=maxlatitude, minlongitude=minlongitude, maxlongitude=maxlongitude,
latitude=latitude, longitude=longitude, minradius=minradius, maxradius=maxradius, level='channel',
includerestricted=includerestricted)
else:
inv = client.get_stations(network=network, station=station, starttime=starttime, endtime=endtime, channel=channel,
minlatitude=minlatitude, maxlatitude=maxlatitude, minlongitude=minlongitude, maxlongitude=maxlongitude,
latitude=latitude, longitude=longitude, minradius=minradius, maxradius=maxradius, level='channel',
includerestricted=includerestricted)
init_flag= True
self.add_stationxml(inv)
try:
self.inv +=inv
except:
self.inv = inv
return
def get_limits_lonlat(self):
"""get the geographical limits of the stations
"""
staLst = self.waveforms.list()
minlat = 90.
maxlat = -90.
minlon = 360.
maxlon = 0.
for staid in staLst:
lat, elv, lon = self.waveforms[staid].coordinates.values()
if lon<0:
lon += 360.
minlat = min(lat, minlat)
maxlat = max(lat, maxlat)
minlon = min(lon, minlon)
maxlon = max(lon, maxlon)
print 'latitude range: ', minlat, '-', maxlat, 'longitude range:', minlon, '-', maxlon
self.minlat = minlat
self.maxlat = maxlat
self.minlon = minlon
self.maxlon = maxlon
return
def _get_basemap(self, projection='lambert', geopolygons=None, resolution='i', blon=0., blat=0.):
"""Get basemap for plotting results
"""
# fig=plt.figure(num=None, figsize=(12, 12), dpi=80, facecolor='w', edgecolor='k')
try:
minlon = self.minlon-blon
maxlon = self.maxlon+blon
minlat = self.minlat-blat
maxlat = self.maxlat+blat
except AttributeError:
self.get_limits_lonlat()
minlon = self.minlon-blon
maxlon = self.maxlon+blon
minlat = self.minlat-blat
maxlat = self.maxlat+blat
lat_centre = (maxlat+minlat)/2.0
lon_centre = (maxlon+minlon)/2.0
if projection == 'merc':
m = Basemap(projection='merc', llcrnrlat=minlat-5., urcrnrlat=maxlat+5., llcrnrlon=minlon-5.,
urcrnrlon=maxlon+5., lat_ts=20, resolution=resolution)
m.drawparallels(np.arange(-80.0,80.0,5.0), labels=[1,0,0,1])
m.drawmeridians(np.arange(-170.0,170.0,5.0), labels=[1,0,0,1])
m.drawstates(color='g', linewidth=2.)
elif projection == 'global':
m = Basemap(projection='ortho',lon_0=lon_centre, lat_0=lat_centre, resolution=resolution)
elif projection == 'regional_ortho':
m1 = Basemap(projection='ortho', lon_0=minlon, lat_0=minlat, resolution='l')
m = Basemap(projection='ortho', lon_0=minlon, lat_0=minlat, resolution=resolution,\
llcrnrx=0., llcrnry=0., urcrnrx=m1.urcrnrx/mapfactor, urcrnry=m1.urcrnry/3.5)
m.drawparallels(np.arange(-80.0,80.0,10.0), labels=[1,0,0,0], linewidth=2, fontsize=20)
m.drawmeridians(np.arange(-170.0,170.0,10.0), linewidth=2)
elif projection=='lambert':
distEW, az, baz = obspy.geodetics.gps2dist_azimuth(minlat, minlon, minlat, maxlon) # distance is in m
distNS, az, baz = obspy.geodetics.gps2dist_azimuth(minlat, minlon, maxlat+2., minlon) # distance is in m
m = Basemap(width = distEW, height=distNS, rsphere=(6378137.00,6356752.3142), resolution='l', projection='lcc',\
lat_1=minlat, lat_2=maxlat, lon_0=lon_centre, lat_0=lat_centre+1)
m.drawparallels(np.arange(-80.0,80.0,10.0), linewidth=1, dashes=[2,2], labels=[1,1,0,0], fontsize=15)
m.drawmeridians(np.arange(-170.0,170.0,10.0), linewidth=1, dashes=[2,2], labels=[0,0,1,0], fontsize=15)
m.drawcoastlines(linewidth=1.0)
m.drawcountries(linewidth=1.)
# m.fillcontinents(lake_color='#99ffff',zorder=0.2)
m.drawmapboundary(fill_color="white")
try:
geopolygons.PlotPolygon(inbasemap=m)
except:
pass
return m
def plot_stations(self, projection='lambert', geopolygons=None, showfig=True, blon=.5, blat=0.5):
"""plot station map
==============================================================================
::: input parameters :::
projection - type of geographical projection
geopolygons - geological polygons for plotting
blon, blat - extending boundaries in longitude/latitude
showfig - show figure or not
==============================================================================
"""
staLst = self.waveforms.list()
stalons = np.array([])
stalats = np.array([])
for staid in staLst:
stla, evz, stlo = self.waveforms[staid].coordinates.values()
stalons = np.append(stalons, stlo); stalats=np.append(stalats, stla)
m = self._get_basemap(projection=projection, geopolygons=geopolygons, blon=blon, blat=blat)
m.etopo()
# m.shadedrelief()
stax, stay = m(stalons, stalats)
m.plot(stax, stay, '^', markersize=5)
# plt.title(str(self.period)+' sec', fontsize=20)
if showfig:
plt.show()
return
def wsac_xcorr(self, netcode1, stacode1, netcode2, stacode2, chan1, chan2, outdir='.', pfx='COR'):
"""Write cross-correlation data from ASDF to sac file
==============================================================================
::: input parameters :::
netcode1, stacode1, chan1 - network/station/channel name for station 1
netcode2, stacode2, chan2 - network/station/channel name for station 2
outdir - output directory
pfx - prefix
::: output :::
e.g. outdir/COR/TA.G12A/COR_TA.G12A_BHT_TA.R21A_BHT.SAC
==============================================================================
"""
subdset = self.auxiliary_data.NoiseXcorr[netcode1][stacode1][netcode2][stacode2][chan1][chan2]
sta1 = self.waveforms[netcode1+'.'+stacode1].StationXML.networks[0].stations[0]
sta2 = self.waveforms[netcode2+'.'+stacode2].StationXML.networks[0].stations[0]
xcorr_sacheader = xcorr_sacheader_default.copy()
xcorr_sacheader['kuser0'] = netcode1
xcorr_sacheader['kevnm'] = stacode1
xcorr_sacheader['knetwk'] = netcode2
xcorr_sacheader['kstnm'] = stacode2
xcorr_sacheader['kcmpnm'] = chan1+chan2
xcorr_sacheader['evla'] = sta1.latitude
xcorr_sacheader['evlo'] = sta1.longitude
xcorr_sacheader['stla'] = sta2.latitude
xcorr_sacheader['stlo'] = sta2.longitude
xcorr_sacheader['dist'] = subdset.parameters['dist']
xcorr_sacheader['az'] = subdset.parameters['az']
xcorr_sacheader['baz'] = subdset.parameters['baz']
xcorr_sacheader['b'] = subdset.parameters['b']
xcorr_sacheader['e'] = subdset.parameters['e']
xcorr_sacheader['delta'] = subdset.parameters['delta']
xcorr_sacheader['npts'] = subdset.parameters['npts']
xcorr_sacheader['user0'] = subdset.parameters['stackday']
sacTr = obspy.io.sac.sactrace.SACTrace(data=subdset.data.value, **xcorr_sacheader)
if not os.path.isdir(outdir+'/'+pfx+'/'+netcode1+'.'+stacode1):
os.makedirs(outdir+'/'+pfx+'/'+netcode1+'.'+stacode1)
sacfname = outdir+'/'+pfx+'/'+netcode1+'.'+stacode1+'/'+ \
pfx+'_'+netcode1+'.'+stacode1+'_'+chan1+'_'+netcode2+'.'+stacode2+'_'+chan2+'.SAC'
sacTr.write(sacfname)
return
def wsac_xcorr_all(self, netcode1, stacode1, netcode2, stacode2, outdir='.', pfx='COR'):
"""Write all components of cross-correlation data from ASDF to sac file
==============================================================================
::: input parameters :::
netcode1, stacode1 - network/station name for station 1
netcode2, stacode2 - network/station name for station 2
outdir - output directory
pfx - prefix
::: output :::
e.g. outdir/COR/TA.G12A/COR_TA.G12A_BHT_TA.R21A_BHT.SAC
==============================================================================
"""
subdset = self.auxiliary_data.NoiseXcorr[netcode1][stacode1][netcode2][stacode2]
channels1 = subdset.list()
channels2 = subdset[channels1[0]].list()
for chan1 in channels1:
for chan2 in channels2:
self.wsac_xcorr(netcode1=netcode1, stacode1=stacode1, netcode2=netcode2,
stacode2=stacode2, chan1=chan1, chan2=chan2, outdir=outdir, pfx=pfx)
return
def get_xcorr_trace(self, netcode1, stacode1, netcode2, stacode2, chan1, chan2):
"""Get one single cross-correlation trace
==============================================================================
::: input parameters :::
netcode1, stacode1, chan1 - network/station/channel name for station 1
netcode2, stacode2, chan2 - network/station/channel name for station 2
::: output :::
obspy trace
==============================================================================
"""
subdset = self.auxiliary_data.NoiseXcorr[netcode1][stacode1][netcode2][stacode2][chan1][chan2]
evla, evz, evlo = self.waveforms[netcode1+'.'+stacode1].coordinates.values()
stla, stz, stlo = self.waveforms[netcode2+'.'+stacode2].coordinates.values()
tr = obspy.core.Trace()
tr.data = subdset.data.value
tr.stats.sac = {}
tr.stats.sac.evla = evla
tr.stats.sac.evlo = evlo
tr.stats.sac.stla = stla
tr.stats.sac.stlo = stlo
tr.stats.sac.kuser0 = netcode1
tr.stats.sac.kevnm = stacode1
tr.stats.network = netcode2
tr.stats.station = stacode2
tr.stats.sac.kcmpnm = chan1+chan2
tr.stats.sac.dist = subdset.parameters['dist']
tr.stats.sac.az = subdset.parameters['az']
tr.stats.sac.baz = subdset.parameters['baz']
tr.stats.sac.b = subdset.parameters['b']
tr.stats.sac.e = subdset.parameters['e']
tr.stats.sac.user0 = subdset.parameters['stackday']
tr.stats.delta = subdset.parameters['delta']
tr.stats.distance = subdset.parameters['dist']*1000.
return tr
def get_xcorr_stream(self, netcode, stacode, chan1, chan2):
st = obspy.Stream()
stalst = self.waveforms.list()
for staid in stalst:
netcode2, stacode2 \
= staid.split('.')
try:
st += self.get_xcorr_trace(netcode1=netcode, stacode1=stacode, netcode2=netcode2, stacode2=stacode2, chan1=chan1, chan2=chan2)
except KeyError:
try:
st += self.get_xcorr_trace(netcode1=netcode2, stacode1=stacode2, netcode2=netcode, stacode2=stacode, chan1=chan2, chan2=chan1)
except KeyError:
pass
return st
def read_xcorr(self, datadir, pfx='COR', fnametype=1, inchannels=None, verbose=True):
"""Read cross-correlation data in ASDF database
===========================================================================================================
::: input parameters :::
datadir - data directory
pfx - prefix
inchannels - input channels, if None, will read channel information from obspy inventory
fnametype - input sac file name type
=1: datadir/2011.JAN/COR/TA.G12A/COR_TA.G12A_BHZ_TA.R21A_BHZ.SAC
=2: datadir/2011.JAN/COR/G12A/COR_G12A_R21A.SAC
=3: datadir/2011.JAN/COR/G12A/COR_G12A_BHZ_R21A_BHZ.SAC
-----------------------------------------------------------------------------------------------------------
::: output :::
ASDF path : self.auxiliary_data.NoiseXcorr[netcode1][stacode1][netcode2][stacode2][chan1][chan2]
===========================================================================================================
"""
staLst = self.waveforms.list()
# main loop for station pairs
if inchannels != None:
try:
if not isinstance(inchannels[0], obspy.core.inventory.channel.Channel):
channels = []
for inchan in inchannels:
channels.append(obspy.core.inventory.channel.Channel(code=inchan, location_code='',
latitude=0, longitude=0, elevation=0, depth=0) )
else:
channels = inchannels
except:
inchannels = None
for staid1 in staLst:
for staid2 in staLst:
netcode1, stacode1 = staid1.split('.')
netcode2, stacode2 = staid2.split('.')
if staid1 >= staid2:
continue
if fnametype==2 and not os.path.isfile(datadir+'/'+pfx+'/'+staid1+'/'+pfx+'_'+staid1+'_'+staid2+'.SAC'):
continue
if inchannels==None:
channels1 = self.waveforms[staid1].StationXML.networks[0].stations[0].channels
channels2 = self.waveforms[staid2].StationXML.networks[0].stations[0].channels
else:
channels1 = channels
channels2 = channels
skipflag = False
for chan1 in channels1:
if skipflag:
break
for chan2 in channels2:
if fnametype == 1:
fname = datadir+'/'+pfx+'/'+staid1+'/'+pfx+'_'+staid1+'_'+chan1.code+'_'\
+staid2+'_'+chan2.code+'.SAC'
elif fnametype == 2:
fname = datadir+'/'+pfx+'/'+stacode1+'/'+pfx+'_'+stacode1+'_'+stacode2+'.SAC'
#----------------------------------------------------------
elif fnametype == 3:
fname = datadir+'/'+pfx+'/'+stacode1+'/'+pfx+'_'+stacode1+'_'+chan1.code+'_'\
+stacode2+'_'+chan2.code+'.SAC'
#----------------------------------------------------------
try:
tr = obspy.core.read(fname)[0]
except IOError:
skipflag= True
break
# write cross-correlation header information
xcorr_header = xcorr_header_default.copy()
xcorr_header['b'] = tr.stats.sac.b
xcorr_header['e'] = tr.stats.sac.e
xcorr_header['netcode1']= netcode1
xcorr_header['netcode2']= netcode2
xcorr_header['stacode1']= stacode1
xcorr_header['stacode2']= stacode2
xcorr_header['npts'] = tr.stats.npts
xcorr_header['delta'] = tr.stats.delta
xcorr_header['stackday']= tr.stats.sac.user0
try:
xcorr_header['dist']= tr.stats.sac.dist
xcorr_header['az'] = tr.stats.sac.az
xcorr_header['baz'] = tr.stats.sac.baz
except AttributeError:
lon1 = self.waveforms[staid1].StationXML.networks[0].stations[0].longitude
lat1 = self.waveforms[staid1].StationXML.networks[0].stations[0].latitude
lon2 = self.waveforms[staid2].StationXML.networks[0].stations[0].longitude
lat2 = self.waveforms[staid2].StationXML.networks[0].stations[0].latitude
dist, az, baz = obspy.geodetics.gps2dist_azimuth(lat1, lon1, lat2, lon2)
dist = dist/1000.
xcorr_header['dist']= dist
xcorr_header['az'] = az
xcorr_header['baz'] = baz
staid_aux = netcode1+'/'+stacode1+'/'+netcode2+'/'+stacode2
xcorr_header['chan1'] = chan1.code
xcorr_header['chan2'] = chan2.code
self.add_auxiliary_data(data=tr.data, data_type='NoiseXcorr', path=staid_aux+'/'+chan1.code+'/'+chan2.code, parameters=xcorr_header)
if verbose and not skipflag:
print 'reading xcorr data: '+netcode1+'.'+stacode1+'_'+netcode2+'.'+stacode2
return
def compute_xcorr(self, datadir, startdate, enddate, chans=['LHZ', 'LHE', 'LHN'], \
fskipxcorr = 0, ftlen = True, tlen = 84000., mintlen = 20000., sps = 1., lagtime = 3000., CorOutflag = 0, \
fprcs = False, fastfft=True, parallel=True, nprocess=None, subsize=1000):
"""
compute ambient noise cross-correlation given preprocessed amplitude and phase files
=================================================================================================================
::: input parameters :::
datadir - directory including data and output
startdate/enddate - start/end date for computation
chans - channel list