-
Notifications
You must be signed in to change notification settings - Fork 29
/
Copy pathMPI_ICESat2_ATL03.py
2449 lines (2376 loc) · 164 KB
/
MPI_ICESat2_ATL03.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
#!/usr/bin/env python
u"""
MPI_ICESat2_ATL03.py (05/2024)
Read ICESat-2 ATL03 and ATL09 data files to calculate average segment surfaces
ATL03 datasets: Global Geolocated Photons
ATL09 datasets: Atmospheric Characteristics
CALLING SEQUENCE:
mpiexec -np 6 python MPI_ICESat2_ATL03.py ATL03_file ATL09_file
COMMAND LINE OPTIONS:
-O X, --output X: Name and path of output file
-V, --verbose: Verbose output to track progress
-M X, --mode X: Permission mode of files created
REQUIRES MPI PROGRAM
MPI: standardized and portable message-passing system
https://www.open-mpi.org/
http://mpitutorial.com/
PYTHON DEPENDENCIES:
numpy: Scientific Computing Tools For Python
https://numpy.org
https://numpy.org/doc/stable/user/numpy-for-matlab-users.html
scipy: Scientific Tools for Python
https://docs.scipy.org/doc/
mpi4py: MPI for Python
http://pythonhosted.org/mpi4py/
http://mpi4py.readthedocs.org/en/stable/
h5py: Python interface for Hierarchal Data Format 5 (HDF5)
https://h5py.org
http://docs.h5py.org/en/stable/mpi.html
scikit-learn: Machine Learning in Python
http://scikit-learn.org/stable/index.html
https://github.com/scikit-learn/scikit-learn
timescale: Python tools for time and astronomical calculations
https://pypi.org/project/timescale/
PROGRAM DEPENDENCIES:
fit.py: Utilities for calculating fits from ATL03 Geolocated Photon Data
utilities.py: download and management utilities for syncing files
classify_photons.py: Yet Another Photon Classifier for Geolocated Photon Data
UPDATE HISTORY:
Updated 05/2024: use wrapper to importlib for optional dependencies
Updated 04/2024: use timescale for temporal operations
Updated 03/2024: use pathlib to define and operate on paths
Updated 12/2022: single implicit import of altimetry tools
Updated 06/2022: update classify photons to match current GSFC version
Updated 05/2022: use argparse descriptions within sphinx documentation
Updated 10/2021: using python logging for handling verbose output
do not use possible TEP photons in photon classification calculation
added parsing for converting file lines to arguments
use reference photon geolocation as default in case of no valid segments
Updated 08/2021: update classify photons to match current GSFC version
Updated 05/2021: add photon classifier based on GSFC YAPC algorithms
move surface fit operations into separate module
Updated 02/2021: replaced numpy bool/int to prevent deprecation warnings
Updated 01/2021: time utilities for converting times from JD and to decimal
Updated 12/2020: H5py deprecation warning change to use make_scale
Updated 10/2020: using argparse to set parameters
Updated 09/2020: using reference photon delta time to interpolate ATL09
Updated 08/2020: using convert delta time function to convert to Julian days
Updated 07/2020: "re-tiding" is no longer unnecessary
Updated 06/2020: verify that complementary beam pair is in list of beams
set masks of output arrays after reading from HDF5
add additional beam check within heights groups
Updated 10/2019: changing Y/N flags to True/False
Updated 09/2019: adding segment quality summary variable
Updated 04/2019: updated backup algorithm for when the surface fit fails
estimate both mean and median first photon bias corrections
estimate both mean and median transmit pulse shape corrections
Updated 03/2019: extract a set of ATL09 parameters for each ATL03 segment_ID
Updated 02/2019: procedures following ATBD for first ATL03 release
Written 05/2017
"""
from __future__ import print_function, division
import sys
import os
import re
import logging
import pathlib
import argparse
import datetime
import warnings
import numpy as np
import scipy.signal
import scipy.interpolate
import icesat2_toolkit as is2tk
import timescale.time
# attempt imports
h5py = is2tk.utilities.import_dependency('h5py')
MPI = is2tk.utilities.import_dependency('mpi4py.MPI')
cluster = is2tk.utilities.import_dependency('sklearn.cluster')
yapc = is2tk.utilities.import_dependency('yapc')
# PURPOSE: keep track of MPI threads
def info(rank, size):
logging.info(f'Rank {rank+1:d} of {size:d}')
logging.info(f'module name: {__name__}')
if hasattr(os, 'getppid'):
logging.info(f'parent process: {os.getppid():d}')
logging.info(f'process id: {os.getpid():d}')
# PURPOSE: create argument parser
def arguments():
parser = argparse.ArgumentParser(
description="""Read ICESat-2 ATL03 and ATL09 data files to calculate
average segment surfaces
""",
fromfile_prefix_chars="@"
)
parser.convert_arg_line_to_args = is2tk.utilities.convert_arg_line_to_args
# command line parameters
# first file listed contains the ATL03 file
# second file listed is the associated ATL09 file
parser.add_argument('ATL03',
type=pathlib.Path, nargs='?',
help='ICESat-2 ATL03 file to run')
parser.add_argument('ATL09',
type=pathlib.Path, nargs='?',
help='ICESat-2 ATL09 file to run')
# use default output file name
parser.add_argument('--output','-O',
type=pathlib.Path,
help='Name and path of output file')
# verbosity settings
# verbose will output information about each output file
parser.add_argument('--verbose','-V',
default=False, action='store_true',
help='Verbose output of run')
# permissions mode of the local files (number in octal)
parser.add_argument('--mode','-M',
type=lambda x: int(x,base=8), default=0o775,
help='Permissions mode of output files')
# return the parser
return parser
# PURPOSE: read ICESat-2 geolocated photon height data (ATL03)
# and backscatter profiles/atmospheric layer characteristics (ATL09)
# Computes average heights over segments
def main():
# start MPI communicator
comm = MPI.COMM_WORLD
# Read the system arguments listed after the program
parser = arguments()
args,_ = parser.parse_known_args()
# create logger
loglevel = logging.INFO if args.verbose else logging.CRITICAL
logging.basicConfig(level=loglevel)
# output module information for process
info(comm.rank,comm.size)
if (comm.rank == 0):
logging.info(f'{str(args.ATL03)} -->')
# compile regular expression operator for extracting data from ATL03 files
rx1 = re.compile(r'(processed_)?(ATL\d+)_(\d{4})(\d{2})(\d{2})(\d{2})(\d{2})'
r'(\d{2})_(\d{4})(\d{2})(\d{2})_(\d{3})_(\d{2})(.*?).h5$')
# universal variables
# speed of light
c = 299792458.0
# associated beam pairs
associated_beam_pair = dict(gt1l='gt1r',gt1r='gt1l',gt2l='gt2r',gt2r='gt2l',
gt3l='gt3r',gt3r='gt3l')
# read ICESat-2 ATL03 HDF5 files (extract base parameters)
SUB,PRD,YY,MM,DD,HH,MN,SS,TRK,CYCL,GRAN,RL,VERS,AUX=rx1.findall(args.ATL03).pop()
# Open the HDF5 file for reading
fileID = h5py.File(args.ATL03, 'r', driver='mpio', comm=comm)
# read each input beam within the file
IS2_atl03_beams = []
for gtx in [k for k in fileID.keys() if bool(re.match(r'gt\d[lr]',k))]:
# check if subsetted beam contains data
# check in both the geolocation and heights groups
try:
fileID[gtx]['geolocation']['segment_id']
fileID[gtx]['heights']['delta_time']
except KeyError:
pass
else:
IS2_atl03_beams.append(gtx)
# number of GPS seconds between the GPS epoch
# and ATLAS Standard Data Product (SDP) epoch
atlas_sdp_gps_epoch = fileID['ancillary_data']['atlas_sdp_gps_epoch'][:]
# which TEP to use for a given spot (convert to 0-based index)
tep_valid_spot = fileID['ancillary_data']['tep']['tep_valid_spot'][:] - 1
tep_pce = ['pce1_spot1','pce2_spot3']
# valid range of times for each TEP histogram
tep_range_prim = fileID['ancillary_data']['tep']['tep_range_prim'][:]
# save tep parameters for a given beam
tep = {}
# variables of interest for generating corrected elevation estimates
Segment_ID = {}
Segment_Index_begin = {}
Segment_PE_count = {}
Segment_Distance = {}
Segment_Length = {}
Segment_Background = {}
# fit parameters
Segment_delta_time = {}
Segment_Height = {}
Segment_Land_Ice = {}
Segment_dH_along = {}
Segment_dH_across = {}
Segment_Height_Error = {}
Segment_Land_Ice_Error = {}
Segment_dH_along_Error = {}
Segment_dH_across_Error = {}
Segment_Mean_Median = {}
Segment_X_atc = {}
Segment_X_spread = {}
Segment_Y_atc = {}
Segment_sigma_geo = {}
Segment_Longitude = {}
Segment_Latitude = {}
Segment_N_Fit = {}
Segment_Window = {}
Segment_RDE = {}
Segment_SNR = {}
Segment_Photon_SNR = {}
Segment_Summary = {}
Segment_Iterations = {}
Segment_Clusters = {}
Segment_Source = {}
Segment_Pulses = {}
# correction parameters
FPB_mean_corr = {}
FPB_mean_sigma = {}
FPB_median_corr = {}
FPB_median_sigma = {}
mean_dead_time = {}
FPB_n_corr = {}
FPB_cal_corr = {}
TPS_mean_corr = {}
TPS_median_corr = {}
# for each input beam within the file
for gtx in sorted(IS2_atl03_beams):
logging.info(gtx) if (comm.rank == 0) else None
# beam type (weak versus strong) for time
atlas_beam_type = fileID[gtx].attrs['atlas_beam_type'].decode('utf-8')
n_pixels = 16.0 if (atlas_beam_type == "strong") else 4.0
# ATL03 Segment ID
Segment_ID[gtx] = fileID[gtx]['geolocation']['segment_id'][:]
# number of valid overlapping ATL03 segments
n_seg = len(Segment_ID[gtx]) - 1
# number of photon events
n_pe, = fileID[gtx]['heights']['delta_time'].shape
# first photon in the segment (convert to 0-based indexing)
Segment_Index_begin[gtx] = fileID[gtx]['geolocation']['ph_index_beg'][:] - 1
# number of photon events in the segment
Segment_PE_count[gtx] = fileID[gtx]['geolocation']['segment_ph_cnt'][:]
# along-track distance for each ATL03 segment
Segment_Distance[gtx] = fileID[gtx]['geolocation']['segment_dist_x'][:]
# along-track length for each ATL03 segment
Segment_Length[gtx] = fileID[gtx]['geolocation']['segment_length'][:]
# ocean tide
fv = fileID[gtx]['geophys_corr']['tide_ocean'].attrs['_FillValue']
tide_ocean = np.ma.array(fileID[gtx]['geophys_corr']['tide_ocean'][:],
fill_value=fv)
tide_ocean.mask = tide_ocean.data == tide_ocean.fill_value
# interpolate background photon rate based on 50-shot summation
background_delta_time = fileID[gtx]['bckgrd_atlas']['delta_time'][:]
SPL = scipy.interpolate.UnivariateSpline(background_delta_time,
fileID[gtx]['bckgrd_atlas']['bckgrd_rate'][:],k=3,s=0)
Segment_Background[gtx] = SPL(fileID[gtx]['geolocation']['delta_time'][:])
# ATLAS spot number for beam in current orientation
spot = int(fileID[gtx].attrs['atlas_spot_number'])
# get ATLAS impulse response variables for the transmitter echo path (TEP)
tep1,tep2 = ('atlas_impulse_response','tep_histogram')
# get appropriate transmitter-echo-path histogram for spot
associated_pce = tep_valid_spot[spot-1]
pce = tep_pce[associated_pce]
# delta time of TEP histogram
tep_tod, = fileID[tep1][pce][tep2]['tep_tod'][:]
# truncate tep to primary histogram (reflection 43-50 ns)
# and extract signal tep from noise tep. calculate width of tep
# ATL03 recommends subsetting between 15-30 ns to avoid secondary
tep_hist_time = np.copy(fileID[tep1][pce][tep2]['tep_hist_time'][:])
tep_hist = np.copy(fileID[tep1][pce][tep2]['tep_hist'][:])
t_TX,p_TX,W_TX,FWHM,TXs,TXe = is2tk.fit.extract_tep_histogram(
tep_hist_time, tep_hist, tep_range_prim)
# save tep information and statistics
tep[gtx] = {}
tep[gtx]['pce'] = pce
tep[gtx]['tep_tod'] = tep_tod
tep[gtx]['tx_start'] = TXs
tep[gtx]['tx_end'] = TXe
tep[gtx]['tx_robust_sprd'] = W_TX
tep[gtx]['sigma_tx'] = FWHM
# channel dead time and first photon bias table for beam
cal1,cal2 = ('ancillary_data','calibrations')
channel_dead_time = fileID[cal1][cal2]['dead_time'][gtx]['dead_time'][:]
mean_dead_time[gtx] = np.mean(channel_dead_time)
fpb_dead_time = fileID[cal1][cal2]['first_photon_bias'][gtx]['dead_time'][:]
fpb_strength = fileID[cal1][cal2]['first_photon_bias'][gtx]['strength'][:]
fpb_width = fileID[cal1][cal2]['first_photon_bias'][gtx]['width'][:]
fpb_corr = fileID[cal1][cal2]['first_photon_bias'][gtx]['ffb_corr'][:]
# calculate first photon bias as a function of strength and width
# for the calculated mean dead time of the beam
ndt,ns,nw = np.shape(fpb_corr)
fpb_corr_dead_time = np.zeros((ns,nw))
for s in range(ns):
for w in range(nw):
SPL = scipy.interpolate.UnivariateSpline(fpb_dead_time/1e9,
fpb_corr[:,s,w],k=3,s=0)
fpb_corr_dead_time[s,w] = SPL(mean_dead_time[gtx])
# bivariate spline for estimating first-photon bias using CAL-19
CAL19 = scipy.interpolate.RectBivariateSpline(fpb_strength[0,:],
fpb_width[0,:]/1e9, fpb_corr_dead_time/1e12, kx=1, ky=1)
# allocate for output segment fit data
fill_value = fileID[gtx]['geolocation']['sigma_h'].attrs['_FillValue']
# delta time of fit photons
# use reference photon delta time as initial guess
# will fill in with segment center delta times
Distributed_delta_time = is2tk.fit.segment_mean(
fileID[gtx]['geolocation']['delta_time'][:])
# segment fit heights
Distributed_Height = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_Height.mask = np.ones((n_seg),dtype=bool)
# land ice height corrected for first photon bias and transmit-pulse shape
Distributed_Land_Ice = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_Land_Ice.mask = np.ones((n_seg),dtype=bool)
# segment fit along-track slopes
Distributed_dH_along = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_dH_along.mask = np.ones((n_seg),dtype=bool)
# segment fit height errors
Distributed_Height_Error = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_Height_Error.mask = np.ones((n_seg),dtype=bool)
# land ice height errors (max of fit or first photon bias uncertainties)
Distributed_Land_Ice_Error = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_Land_Ice_Error.mask = np.ones((n_seg),dtype=bool)
# segment fit along-track slope errors
Distributed_dH_along_Error = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_dH_along_Error.mask = np.ones((n_seg),dtype=bool)
# difference between the mean and median of the residuals from fit height
Distributed_Mean_Median = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_Mean_Median.mask = np.ones((n_seg),dtype=bool)
# along-track X coordinates of segment fit
Distributed_X_atc = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_X_atc.mask = np.ones((n_seg),dtype=bool)
# along-track X coordinate spread of points used in segment fit
Distributed_X_spread = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_X_spread.mask = np.ones((n_seg),dtype=bool)
# along-track Y coordinates of segment fit
Distributed_Y_atc = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_Y_atc.mask = np.ones((n_seg),dtype=bool)
# longitude of fit photons
Distributed_Longitude = is2tk.fit.segment_mean(
fileID[gtx]['geolocation']['reference_photon_lon'][:])
# latitude of fit photons
Distributed_Latitude = is2tk.fit.segment_mean(
fileID[gtx]['geolocation']['reference_photon_lat'][:])
# number of photons in fit
Distributed_N_Fit = np.ma.zeros((n_seg),fill_value=-1,dtype=int)
Distributed_N_Fit.mask = np.ones((n_seg),dtype=bool)
# size of the window used in the fit
Distributed_Window = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_Window.mask = np.ones((n_seg),dtype=bool)
# robust dispersion estimator
Distributed_RDE = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_RDE.mask = np.ones((n_seg),dtype=bool)
# signal-to-noise ratio
Distributed_SNR = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_SNR.mask = np.ones((n_seg),dtype=bool)
# maximum signal-to-noise ratio from photon classifier
Distributed_Photon_SNR = np.ma.zeros((n_seg),fill_value=0,dtype=int)
Distributed_Photon_SNR.mask = np.ones((n_seg),dtype=bool)
# segment quality summary
Distributed_Summary = np.ma.zeros((n_seg),fill_value=-1,dtype=int)
Distributed_Summary.mask = np.ones((n_seg),dtype=bool)
# number of iterations for fit
Distributed_Iterations = np.ma.zeros((n_seg),fill_value=-1,dtype=int)
Distributed_Iterations.mask = np.ones((n_seg),dtype=bool)
# number of estimated clusters of data
Distributed_Clusters = np.ma.zeros((n_seg),fill_value=0,dtype=int)
Distributed_Clusters.mask = np.ones((n_seg),dtype=bool)
# signal source selection
Distributed_Source = np.ma.zeros((n_seg),fill_value=4,dtype=int)
Distributed_Source.mask = np.ones((n_seg),dtype=bool)
# number of pulses in segment
Distributed_Pulses = np.ma.zeros((n_seg),fill_value=-1,dtype=int)
Distributed_Pulses.mask = np.ones((n_seg),dtype=bool)
# first photon bias estimates
Distributed_FPB_mean_corr = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_FPB_mean_corr.mask = np.ones((n_seg),dtype=bool)
Distributed_FPB_mean_sigma = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_FPB_mean_sigma.mask = np.ones((n_seg),dtype=bool)
Distributed_FPB_median_corr = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_FPB_median_corr.mask = np.ones((n_seg),dtype=bool)
Distributed_FPB_median_sigma = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_FPB_median_sigma.mask = np.ones((n_seg),dtype=bool)
Distributed_FPB_n_corr = np.ma.zeros((n_seg),fill_value=-1,dtype=int)
Distributed_FPB_n_corr.mask = np.ones((n_seg),dtype=bool)
Distributed_FPB_cal_corr = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_FPB_cal_corr.mask = np.ones((n_seg),dtype=bool)
# transmit pulse shape bias estimates
Distributed_TPS_mean_corr = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_TPS_mean_corr.mask = np.ones((n_seg),dtype=bool)
Distributed_TPS_median_corr = np.ma.zeros((n_seg),fill_value=fill_value)
Distributed_TPS_median_corr.mask = np.ones((n_seg),dtype=bool)
# along-track and across-track distance for photon events
x_atc = fileID[gtx]['heights']['dist_ph_along'][:].copy()
y_atc = fileID[gtx]['heights']['dist_ph_across'][:].copy()
# photon event heights
h_ph = fileID[gtx]['heights']['h_ph'][:].copy()
# digital elevation model interpolated to photon events
dem_h = np.zeros((n_pe))
# for each 20m segment
for j,_ in enumerate(Segment_ID[gtx]):
# index for 20m segment j
idx = Segment_Index_begin[gtx][j]
# skip segments with no photon events
if (idx < 0):
continue
# number of photons in 20m segment
cnt = Segment_PE_count[gtx][j]
# add segment distance to along-track coordinates
x_atc[idx:idx+cnt] += Segment_Distance[gtx][j]
# interpolate digital elevation model to photon events
dem_h[idx:idx+cnt] = fileID[gtx]['geophys_corr']['dem_h'][j]
# iterate over ATLAS major frames
photon_mframes = fileID[gtx]['heights']['pce_mframe_cnt'][:].copy()
# background ATLAS group variables are based upon 50-shot summations
# PCE Major Frames are based upon 200-shot summations
pce_mframe_cnt = fileID[gtx]['bckgrd_atlas']['pce_mframe_cnt'][:].copy()
# find unique major frames and their indices within background ATLAS group
# (there will 4 background ATLAS time steps for nearly every major frame)
unique_major_frames,unique_index = np.unique(pce_mframe_cnt,return_index=True)
# number of unique major frames in granule for beam
major_frame_count = len(unique_major_frames)
# height of each telemetry band for a major frame
tlm_height = {}
tlm_height['band1'] = fileID[gtx]['bckgrd_atlas']['tlm_height_band1'][:].copy()
tlm_height['band2'] = fileID[gtx]['bckgrd_atlas']['tlm_height_band2'][:].copy()
# elevation above ellipsoid of each telemetry band for a major frame
tlm_top = {}
tlm_top['band1'] = fileID[gtx]['bckgrd_atlas']['tlm_top_band1'][:].copy()
tlm_top['band2'] = fileID[gtx]['bckgrd_atlas']['tlm_top_band2'][:].copy()
# buffer to telemetry band to set as valid
tlm_buffer = 100.0
# flag denoting photon events as possible TEP
if (int(RL) < 4):
isTEP = np.any((fileID[gtx]['heights']['signal_conf_ph'][:]==-2),axis=1)
else:
isTEP = (fileID[gtx]['heights']['quality_ph'][:] == 3)
# photon event weights
Distributed_Weights = np.zeros((n_pe),dtype=np.float64)
# run for each major frame (distributed over comm.size # of processes)
for iteration in range(comm.rank, major_frame_count, comm.size):
# background atlas index for iteration
idx = unique_index[iteration]
# photon indices for major frame (buffered by 1 frame on each side)
# do not use possible TEP photons in photon classification
i1, = np.nonzero((photon_mframes >= unique_major_frames[iteration]-1) &
(photon_mframes <= unique_major_frames[iteration]+1) &
np.logical_not(isTEP))
# indices for the major frame within the buffered window
i2, = np.nonzero(photon_mframes[i1] == unique_major_frames[iteration])
# sum of telemetry band widths for major frame
h_win_width = 0.0
# check that each telemetry band is close to DEM
for b in ['band1','band2']:
# bottom of the telemetry band for major frame
tlm_bot_band = tlm_top[b][idx] - tlm_height[b][idx]
if np.any((dem_h[i1[i2]] >= (tlm_bot_band-tlm_buffer)) &
(dem_h[i1[i2]] <= (tlm_top[b][idx]+tlm_buffer))):
# add telemetry height to window width
h_win_width += tlm_height[b][idx]
# calculate photon event weights
Distributed_Weights[i1[i2]] = yapc.classify_photons(x_atc[i1],
h_ph[i1], h_win_width, i2, K=0, min_knn=5, min_ph=3,
min_xspread=1.0, min_hspread=0.01, win_x=15.0, win_h=6.0,
method='linear')
# photon event weights
pe_weights = np.zeros((n_pe),dtype=np.float64)
comm.Allreduce(sendbuf=[Distributed_Weights, MPI.DOUBLE], \
recvbuf=[pe_weights, MPI.DOUBLE], op=MPI.SUM)
Distributed_Weights = None
# wait for all distributed processes to finish for beam
comm.Barrier()
# photon event weights scaled to a single byte
weight_ph = np.array(255*pe_weights,dtype=np.uint8)
# verify photon event weights
np.clip(weight_ph, 0, 255, out=weight_ph)
# iterate over valid ATL03 segments
# in ATL03 1-based indexing: invalid == 0
# here in 0-based indexing: invalid == -1
segment_indices, = np.nonzero((Segment_Index_begin[gtx][:-1] >= 0) &
(Segment_Index_begin[gtx][1:] >= 0))
iteration_count = len(segment_indices)
# run for each geoseg (distributed over comm.size # of processes)
for iteration in range(comm.rank, iteration_count, comm.size):
# indice for iteration (can run through a subset of segments)
j = segment_indices[iteration]
# iterate over valid ATL03 segments
# in ATL03 1-based indexing: invalid == 0
# here in 0-based indexing: invalid == -1
if (Segment_Index_begin[gtx][j] >= 0):
# index for segment j
idx = Segment_Index_begin[gtx][j]
# number of photons in segment (use 2 ATL03 segments)
c1 = int(Segment_PE_count[gtx][j])
c2 = int(Segment_PE_count[gtx][j+1])
cnt = c1 + c2
# time of each Photon event (PE)
segment_times = np.copy(fileID[gtx]['heights']['delta_time'][idx:idx+cnt])
# Photon event lat/lon and elevation (re-tided WGS84)
segment_heights = np.copy(h_ph[idx:idx+cnt])
# ATL03 pe heights no longer apply the ocean tide
# and so "re-tiding" is no longer unnecessary
# segment_heights[:c1] += tide_ocean[j]
# segment_heights[c1:] += tide_ocean[j+1]
segment_lats = np.copy(fileID[gtx]['heights']['lat_ph'][idx:idx+cnt])
segment_lons = np.copy(fileID[gtx]['heights']['lon_ph'][idx:idx+cnt])
# Photon event channel and identification
ID_channel = np.copy(fileID[gtx]['heights']['ph_id_channel'][idx:idx+cnt])
ID_pulse = np.copy(fileID[gtx]['heights']['ph_id_pulse'][idx:idx+cnt])
n_pulses = np.unique(ID_pulse).__len__()
frame_number = np.copy(fileID[gtx]['heights']['pce_mframe_cnt'][idx:idx+cnt])
# vertical noise-photon density
background_density = 2.0*n_pulses*Segment_Background[gtx][j]/c
# along-track X and Y coordinates
distance_along_X = np.copy(x_atc[idx:idx+cnt])
distance_along_Y = np.copy(y_atc[idx:idx+cnt])
# check the spread of photons along-track (must be > 20m)
along_X_spread = distance_along_X.max() - distance_along_X.min()
# check confidence level associated with each photon event
# -2: TEP
# -1: Events not associated with a specific surface type
# 0: noise
# 1: buffer but algorithm classifies as background
# 2: low
# 3: medium
# 4: high
# Surface types for signal classification confidence
# 0=Land; 1=Ocean; 2=SeaIce; 3=LandIce; 4=InlandWater
ice_sig_conf = np.copy(fileID[gtx]['heights']['signal_conf_ph'][idx:idx+cnt,3])
ice_sig_low_count = np.count_nonzero(ice_sig_conf > 1)
# indices of TEP classified photons
ice_sig_tep_pe, = np.nonzero(ice_sig_conf == -2)
# photon event weights from photon classifier
segment_weights = weight_ph[idx:idx+cnt]
snr_norm = np.max(segment_weights)
# photon event signal-to-noise ratio from photon classifier
photon_snr = np.array(100.0*segment_weights/snr_norm,dtype=int)
Distributed_Photon_SNR.data[j] = np.copy(snr_norm)
Distributed_Photon_SNR.mask[j] = (snr_norm > 0)
# photon confidence levels from classifier
pe_sig_conf = np.zeros((cnt),dtype=int)
# calculate confidence levels from photon classifier
pe_sig_conf[photon_snr >= 25] = 2
pe_sig_conf[photon_snr >= 60] = 3
pe_sig_conf[photon_snr >= 80] = 4
# copy classification for TEP photons
pe_sig_conf[ice_sig_tep_pe] = -2
pe_sig_low_count = np.count_nonzero(pe_sig_conf > 1)
# check if segment has photon events classified for land ice
# that are at or above low-confidence threshold
# and that the spread of photons is greater than 20m
if (pe_sig_low_count > 10) & (along_X_spread > 20):
# use density-based spatial clustering in segment
db = cluster.DBSCAN(eps=0.5).fit(
np.c_[distance_along_X, segment_heights],
sample_weight=photon_snr)
labels = db.labels_
# number of noise photons
noise_photons = list(labels).count(-1)
noise_cluster = 1 if noise_photons else 0
# number of photon event clusters in segment
n_clusters = len(set(labels)) - noise_cluster
Distributed_Clusters.data[j] = n_clusters
Distributed_Clusters.mask[j] = (n_clusters > 0)
# perform a surface fit procedure
Segment_X = Segment_Distance[gtx][j] + Segment_Length[gtx][j]
valid,fit,centroid = is2tk.fit.try_surface_fit(
distance_along_X, distance_along_Y, segment_heights,
pe_sig_conf, Segment_X, SURF_TYPE='linear', ITERATE=20,
CONFIDENCE=[1,0])
# indices of points used in final iterated fit
ifit = fit['indices'] if valid else None
if bool(valid) & (np.abs(fit['error'][0]) < 20):
Distributed_Height.data[j] = fit['beta'][0]
Distributed_Height.mask[j] = False
Distributed_dH_along.data[j] = fit['beta'][1]
Distributed_dH_along.mask[j] = False
Distributed_Height_Error.data[j] = fit['error'][0]
Distributed_Height_Error.mask[j] = False
Distributed_dH_along_Error.data[j] = fit['error'][1]
Distributed_dH_along_Error.mask[j] = False
# along-track and cross-track coordinates
Distributed_X_atc.data[j] = np.copy(centroid['x'])
Distributed_X_atc.mask[j] = False
Distributed_X_spread.data[j] = np.copy(along_X_spread)
Distributed_X_spread.mask[j] = False
Distributed_Y_atc.data[j] = np.copy(centroid['y'])
Distributed_Y_atc.mask[j] = False
# fit geolocation to the along-track distance of segment
Distributed_delta_time[j] = \
is2tk.fit.fit_geolocation(segment_times[ifit],
distance_along_X[ifit], Distributed_X_atc[j])
Distributed_Longitude[j] = \
is2tk.fit.fit_geolocation(segment_lons[ifit],
distance_along_X[ifit], Distributed_X_atc[j])
Distributed_Latitude[j] = \
is2tk.fit.fit_geolocation(segment_lats[ifit],
distance_along_X[ifit], Distributed_X_atc[j])
# number of photons used in fit
Distributed_N_Fit.data[j] = len(ifit)
Distributed_N_Fit.mask[j] = False
# size of the final window
Distributed_Window.data[j] = np.copy(fit['window'])
Distributed_Window.mask[j] = False
# robust dispersion estimator
Distributed_RDE.data[j] = np.copy(fit['RDE'])
Distributed_RDE.mask[j] = False
# signal to noise ratio
N_BG = background_density*Distributed_Window.data[j]
Distributed_SNR.data[j] = Distributed_N_Fit.data[j]/N_BG
Distributed_SNR.mask[j] = False
# number of iterations used in fit
Distributed_Iterations.data[j] = np.copy(fit['iterations'])
Distributed_Iterations.mask[j] = False
Distributed_Source.data[j] = np.copy(valid)
Distributed_Source.mask[j] = False
Distributed_Pulses.data[j] = np.copy(n_pulses)
Distributed_Pulses.mask[j] = False
# calculate residuals off of fit surface for all data
x_slope = Distributed_dH_along[j]*(distance_along_X-Distributed_X_atc[j])
height_residuals = segment_heights-Distributed_Height[j]-x_slope
temporal_residuals = -2.0*height_residuals/c
# calculate difference between the mean and the median from the fit
Distributed_Mean_Median.data[j] = np.mean(height_residuals[ifit]) - \
np.median(height_residuals[ifit])
Distributed_Mean_Median.mask[j] = False
# calculate flags for quality summary
VPD = Distributed_N_Fit.data[j]/Distributed_Window.data[j]
Distributed_Summary.data[j] = int(
(Distributed_RDE.data[j] >= 1) |
(Distributed_Height_Error.data[j] >= 1) |
(VPD <= (n_pixels/4.0)))
Distributed_Summary.mask[j] = False
# estimate first photon bias corrections
# step-size for histograms (50 ps ~ 7.5mm height)
ii, = np.nonzero((height_residuals >= -Distributed_Window.data[j]) &
(height_residuals <= Distributed_Window.data[j]))
try:
FPB = is2tk.fit.calc_first_photon_bias(
temporal_residuals[ii], n_pulses, n_pixels,
mean_dead_time[gtx], 5e-11, ITERATE=20)
except:
pass
else:
Distributed_FPB_mean_corr.data[j] = -0.5*FPB['mean']*c
Distributed_FPB_mean_corr.mask[j] = False
Distributed_FPB_mean_sigma.data[j] = 0.5*FPB['mean_sigma']*c
Distributed_FPB_mean_sigma.mask[j] = False
Distributed_FPB_median_corr.data[j] = -0.5*FPB['median']*c
Distributed_FPB_median_corr.mask[j] = False
Distributed_FPB_median_sigma.data[j] = 0.5*FPB['median_sigma']*c
Distributed_FPB_median_sigma.mask[j] = False
Distributed_FPB_n_corr.data[j] = np.copy(FPB['count'])
Distributed_FPB_n_corr.mask[j] = False
# first photon bias correction from CAL-19
FPB_calibrated = CAL19.ev(FPB['strength'],FPB['width'])
Distributed_FPB_cal_corr.data[j] = -0.5*FPB_calibrated*c
Distributed_FPB_cal_corr.mask[j] = False
# estimate transmit pulse shape correction
try:
W_RX = 2.0*Distributed_RDE.data[j]/c
dt_W = 2.0*Distributed_Window.data[j]/c
TPS = is2tk.fit.calc_transmit_pulse_shape(t_TX,
p_TX, W_TX, W_RX, dt_W, Distributed_SNR.data[j],
ITERATE=50)
except:
pass
else:
Distributed_TPS_mean_corr.data[j] = 0.5*TPS['mean']*c
Distributed_TPS_mean_corr.mask[j] = False
Distributed_TPS_median_corr.data[j] = 0.5*TPS['median']*c
Distributed_TPS_median_corr.mask[j] = False
# some ATL03 segments will not result in a valid fit
# backup algorithm uses 4 segments to find a valid surface
if (j not in (0,n_seg-2,n_seg-1)) & Distributed_Height.mask[j] & \
(Segment_Index_begin[gtx][j-1] > 0):
# index for segment j
idx = Segment_Index_begin[gtx][j-1]
# number of photons in segment (use 4 ATL03 segments)
c1 = Segment_PE_count[gtx][j-1].astype(int)
c2 = Segment_PE_count[gtx][j].astype(int)
c3 = Segment_PE_count[gtx][j+1].astype(int)
c4 = Segment_PE_count[gtx][j+2].astype(int)
cnt = c1 + c2 + c3 + c4
# time of each Photon event (PE)
segment_times = np.copy(fileID[gtx]['heights']['delta_time'][idx:idx+cnt])
# Photon event lat/lon and elevation (re-tided WGS84)
segment_heights = np.copy(h_ph[idx:idx+cnt])
# ATL03 pe heights no longer apply the ocean tide
# and so "re-tiding" is no longer unnecessary
# segment_heights[:c1] += tide_ocean[j-1]
# segment_heights[c1:c1+c2] += tide_ocean[j]
# segment_heights[c1+c2:c1+c2+c3] += tide_ocean[j+1]
# segment_heights[c1+c2+c3:] += tide_ocean[j+2]
segment_lats = np.copy(fileID[gtx]['heights']['lat_ph'][idx:idx+cnt])
segment_lons = np.copy(fileID[gtx]['heights']['lon_ph'][idx:idx+cnt])
# Photon event channel and identification
ID_channel = np.copy(fileID[gtx]['heights']['ph_id_channel'][idx:idx+cnt])
ID_pulse = np.copy(fileID[gtx]['heights']['ph_id_pulse'][idx:idx+cnt])
n_pulses = np.unique(ID_pulse).__len__()
frame_number = np.copy(fileID[gtx]['heights']['pce_mframe_cnt'][idx:idx+cnt])
# vertical noise-photon density
background_density = 2.0*n_pulses*Segment_Background[gtx][j]/c
# along-track X and Y coordinates
distance_along_X = np.copy(x_atc[idx:idx+cnt])
distance_along_Y = np.copy(y_atc[idx:idx+cnt])
# check the spread of photons along-track (must be > 40m)
along_X_spread = distance_along_X.max() - distance_along_X.min()
# check confidence level associated with each photon event
# -2: TEP
# -1: Events not associated with a specific surface type
# 0: noise
# 1: buffer but algorithm classifies as background
# 2: low
# 3: medium
# 4: high
# Surface types for signal classification confidence
# 0=Land; 1=Ocean; 2=SeaIce; 3=LandIce; 4=InlandWater
ice_sig_conf = np.copy(fileID[gtx]['heights']['signal_conf_ph'][idx:idx+cnt,3])
ice_sig_low_count = np.count_nonzero(ice_sig_conf > 1)
# indices of TEP classified photons
ice_sig_tep_pe, = np.nonzero(ice_sig_conf == -2)
# photon event weights from photon classifier
segment_weights = pe_weights[idx:idx+cnt]
snr_norm = np.max(segment_weights)
# photon event signal-to-noise ratio from photon classifier
photon_snr = np.zeros((cnt),dtype=int)
if (snr_norm > 0):
photon_snr[:] = 100.0*segment_weights/snr_norm
# copy signal to noise ratio for segment
Distributed_Photon_SNR.data[j] = np.copy(snr_norm)
Distributed_Photon_SNR.mask[j] = (snr_norm > 0)
# photon confidence levels from classifier
pe_sig_conf = np.zeros((cnt),dtype=int)
# calculate confidence levels from photon classifier
pe_sig_conf[photon_snr >= 25] = 2
pe_sig_conf[photon_snr >= 60] = 3
pe_sig_conf[photon_snr >= 80] = 4
# copy classification for TEP photons
pe_sig_conf[ice_sig_tep_pe] = -2
pe_sig_low_count = np.count_nonzero(pe_sig_conf > 1)
# check if segment has photon events classified for land ice
# that are at or above low-confidence threshold
# and that the spread of photons is greater than 40m
if (pe_sig_low_count > 10) & (along_X_spread > 40):
# use density-based spatial clustering in segment
db = cluster.DBSCAN(eps=0.5).fit(
np.c_[distance_along_X, segment_heights],
sample_weight=photon_snr)
labels = db.labels_
# number of noise photons
noise_photons = list(labels).count(-1)
noise_cluster = 1 if noise_photons else 0
# number of photon event clusters in segment
n_clusters = len(set(labels)) - noise_cluster
Distributed_Clusters.data[j] = n_clusters
Distributed_Clusters.mask[j] = (n_clusters > 0)
# perform a surface fit procedure
Segment_X = Segment_Distance[gtx][j] + Segment_Length[gtx][j]
valid,fit,centroid = is2tk.fit.try_surface_fit(
distance_along_X, distance_along_Y, segment_heights,
pe_sig_conf, Segment_X, SURF_TYPE='quadratic',
ITERATE=20, CONFIDENCE=[0])
# indices of points used in final iterated fit
ifit = fit['indices'] if valid else None
if bool(valid) & (np.abs(fit['error'][0]) < 20):
Distributed_Height.data[j] = fit['beta'][0]
Distributed_Height.mask[j] = False
Distributed_dH_along.data[j] = fit['beta'][1]
Distributed_dH_along.mask[j] = False
Distributed_Height_Error.data[j] = fit['error'][0]
Distributed_Height_Error.mask[j] = False
Distributed_dH_along_Error.data[j] = fit['error'][1]
Distributed_dH_along_Error.mask[j] = False
# along-track and cross-track coordinates
Distributed_X_atc.data[j] = np.copy(centroid['x'])
Distributed_X_atc.mask[j] = False
Distributed_X_spread.data[j] = np.copy(along_X_spread)
Distributed_X_spread.mask[j] = False
Distributed_Y_atc.data[j] = np.copy(centroid['y'])
Distributed_Y_atc.mask[j] = False
# fit geolocation to the along-track distance of segment
Distributed_delta_time[j] = \
is2tk.fit.fit_geolocation(segment_times[ifit],
distance_along_X[ifit], Distributed_X_atc[j])
Distributed_Longitude[j] = \
is2tk.fit.fit_geolocation(segment_lons[ifit],
distance_along_X[ifit], Distributed_X_atc[j])
Distributed_Latitude[j] = \
is2tk.fit.fit_geolocation(segment_lats[ifit],
distance_along_X[ifit], Distributed_X_atc[j])
# number of photons used in fit
Distributed_N_Fit.data[j] = len(ifit)
Distributed_N_Fit.mask[j] = False
# size of the final window
Distributed_Window.data[j] = np.copy(fit['window'])
Distributed_Window.mask[j] = False
# robust dispersion estimator
Distributed_RDE.data[j] = np.copy(fit['RDE'])
Distributed_RDE.mask[j] = False
# signal to noise ratio
N_BG = background_density*Distributed_Window.data[j]
Distributed_SNR.data[j] = Distributed_N_Fit.data[j]/N_BG
Distributed_SNR.mask[j] = False
# number of iterations used in fit
Distributed_Iterations.data[j] = np.copy(fit['iterations'])
Distributed_Iterations.mask[j] = False
Distributed_Source.data[j] = 2 + np.copy(valid)
Distributed_Source.mask[j] = False
Distributed_Pulses.data[j] = np.copy(n_pulses)
Distributed_Pulses.mask[j] = False
# calculate residuals off of fit surface for all data
x_slope = Distributed_dH_along[j]*(distance_along_X-Distributed_X_atc[j])
height_residuals = segment_heights-Distributed_Height[j]-x_slope
temporal_residuals = -2.0*height_residuals/c
# calculate difference between the mean and the median from the fit
Distributed_Mean_Median.data[j] = np.mean(height_residuals[ifit]) - \
np.median(height_residuals[ifit])
Distributed_Mean_Median.mask[j] = False
# calculate flags for quality summary
VPD = Distributed_N_Fit.data[j]/Distributed_Window.data[j]
Distributed_Summary.data[j] = int(
(Distributed_RDE.data[j] >= 1) |
(Distributed_Height_Error.data[j] >= 1) |
(VPD <= (n_pixels/4.0)))
Distributed_Summary.mask[j] = False
# estimate first photon bias corrections
# step-size for histograms (50 ps ~ 7.5mm height)
try:
ii, = np.nonzero((height_residuals >= -Distributed_Window.data[j]) &
(height_residuals <= Distributed_Window.data[j]))
FPB = is2tk.fit.calc_first_photon_bias(
temporal_residuals[ii], n_pulses, n_pixels,
mean_dead_time[gtx], 5e-11, ITERATE=20)
except:
pass
else:
Distributed_FPB_mean_corr.data[j] = -0.5*FPB['mean']*c
Distributed_FPB_mean_corr.mask[j] = False
Distributed_FPB_mean_sigma.data[j] = 0.5*FPB['mean_sigma']*c
Distributed_FPB_mean_sigma.mask[j] = False
Distributed_FPB_median_corr.data[j] = -0.5*FPB['median']*c
Distributed_FPB_median_corr.mask[j] = False
Distributed_FPB_median_sigma.data[j] = 0.5*FPB['median_sigma']*c
Distributed_FPB_median_sigma.mask[j] = False
Distributed_FPB_n_corr.data[j] = np.copy(FPB['count'])
Distributed_FPB_n_corr.mask[j] = False
# first photon bias correction from CAL-19
FPB_calibrated = CAL19.ev(FPB['strength'],FPB['width'])
Distributed_FPB_cal_corr.data[j] = -0.5*FPB_calibrated*c
Distributed_FPB_cal_corr.mask[j] = False
# estimate transmit pulse shape correction
try:
W_RX = 2.0*Distributed_RDE.data[j]/c
dt_W = 2.0*Distributed_Window.data[j]/c
TPS = is2tk.fit.calc_transmit_pulse_shape(t_TX,
p_TX, W_TX, W_RX, dt_W, Distributed_SNR.data[j],
ITERATE=50)
except:
pass
else:
Distributed_TPS_mean_corr.data[j] = 0.5*TPS['mean']*c
Distributed_TPS_mean_corr.mask[j] = False
Distributed_TPS_median_corr.data[j] = 0.5*TPS['median']*c
Distributed_TPS_median_corr.mask[j] = False
# if there is a valid land ice height
if (~Distributed_Height.mask[j]):
# land ice height corrected for first photon bias and transmit-pulse shape
# segment heights have already been "re-tided"
Distributed_Land_Ice.data[j] = Distributed_Height.data[j] + \
Distributed_FPB_median_corr.data[j] + Distributed_TPS_median_corr.data[j]
Distributed_Land_Ice.mask[j] = False
# land ice height errors (max of fit or first photon bias uncertainties)
Distributed_Land_Ice_Error.data[j] = np.sqrt(np.max([
Distributed_Height_Error.data[j]**2,
Distributed_FPB_median_sigma.data[j]**2]))
Distributed_Land_Ice_Error.mask[j] = False
# communicate output MPI matrices between ranks
# operations are element summations and logical "and" across elements
# delta time of fit photons
Segment_delta_time[gtx] = np.zeros((n_seg),fill_value=fill_value)
comm.Allreduce(sendbuf=[Distributed_delta_time, MPI.DOUBLE], \
recvbuf=[Segment_delta_time[gtx], MPI.DOUBLE], op=MPI.SUM)
Distributed_delta_time = None
# segment fit heights
Segment_Height[gtx] = np.ma.zeros((n_seg),fill_value=fill_value)
Segment_Height[gtx].mask = np.ones((n_seg),dtype=bool)
comm.Allreduce(sendbuf=[Distributed_Height.data, MPI.DOUBLE], \
recvbuf=[Segment_Height[gtx].data, MPI.DOUBLE], op=MPI.SUM)
comm.Allreduce(sendbuf=[Distributed_Height.mask, MPI.BOOL], \
recvbuf=[Segment_Height[gtx].mask, MPI.BOOL], op=MPI.LAND)
Distributed_Height = None
# land ice height corrected for first photon bias and transmit-pulse shape
Segment_Land_Ice[gtx] = np.ma.zeros((n_seg),fill_value=fill_value)
Segment_Land_Ice[gtx].mask = np.ones((n_seg),dtype=bool)
comm.Allreduce(sendbuf=[Distributed_Land_Ice.data, MPI.DOUBLE], \
recvbuf=[Segment_Land_Ice[gtx].data, MPI.DOUBLE], op=MPI.SUM)
comm.Allreduce(sendbuf=[Distributed_Land_Ice.mask, MPI.BOOL], \
recvbuf=[Segment_Land_Ice[gtx].mask, MPI.BOOL], op=MPI.LAND)
Distributed_Land_Ice = None
# segment fit along-track slopes
Segment_dH_along[gtx] = np.ma.zeros((n_seg),fill_value=fill_value)
Segment_dH_along[gtx].mask = np.ones((n_seg),dtype=bool)
comm.Allreduce(sendbuf=[Distributed_dH_along.data, MPI.DOUBLE], \
recvbuf=[Segment_dH_along[gtx].data, MPI.DOUBLE], op=MPI.SUM)
comm.Allreduce(sendbuf=[Distributed_dH_along.mask, MPI.BOOL], \
recvbuf=[Segment_dH_along[gtx].mask, MPI.BOOL], op=MPI.LAND)
Distributed_dH_along = None
# segment fit height errors
Segment_Height_Error[gtx] = np.ma.zeros((n_seg),fill_value=fill_value)
Segment_Height_Error[gtx].mask = np.ones((n_seg),dtype=bool)
comm.Allreduce(sendbuf=[Distributed_Height_Error.data, MPI.DOUBLE], \
recvbuf=[Segment_Height_Error[gtx].data, MPI.DOUBLE], op=MPI.SUM)
comm.Allreduce(sendbuf=[Distributed_Height_Error.mask, MPI.BOOL], \
recvbuf=[Segment_Height_Error[gtx].mask, MPI.BOOL], op=MPI.LAND)
Distributed_Height_Error = None
# land ice height errors (max of fit or first photon bias uncertainties)
Segment_Land_Ice_Error[gtx] = np.ma.zeros((n_seg),fill_value=fill_value)
Segment_Land_Ice_Error[gtx].mask = np.ones((n_seg),dtype=bool)
comm.Allreduce(sendbuf=[Distributed_Land_Ice_Error.data, MPI.DOUBLE], \
recvbuf=[Segment_Land_Ice_Error[gtx].data, MPI.DOUBLE], op=MPI.SUM)
comm.Allreduce(sendbuf=[Distributed_Land_Ice_Error.mask, MPI.BOOL], \
recvbuf=[Segment_Land_Ice_Error[gtx].mask, MPI.BOOL], op=MPI.LAND)
Distributed_Land_Ice_Error = None
# segment fit along-track slope errors
Segment_dH_along_Error[gtx] = np.ma.zeros((n_seg),fill_value=fill_value)
Segment_dH_along_Error[gtx].mask = np.ones((n_seg),dtype=bool)
comm.Allreduce(sendbuf=[Distributed_dH_along_Error.data, MPI.DOUBLE], \
recvbuf=[Segment_dH_along_Error[gtx].data, MPI.DOUBLE], op=MPI.SUM)
comm.Allreduce(sendbuf=[Distributed_dH_along_Error.mask, MPI.BOOL], \
recvbuf=[Segment_dH_along_Error[gtx].mask, MPI.BOOL], op=MPI.LAND)
Distributed_dH_along_Error = None
# difference between the mean and median of the residuals from fit height
Segment_Mean_Median[gtx] = np.ma.zeros((n_seg),fill_value=fill_value)
Segment_Mean_Median[gtx].mask = np.ones((n_seg),dtype=bool)
comm.Allreduce(sendbuf=[Distributed_Mean_Median.data, MPI.DOUBLE], \
recvbuf=[Segment_Mean_Median[gtx].data, MPI.DOUBLE], op=MPI.SUM)
comm.Allreduce(sendbuf=[Distributed_Mean_Median.mask, MPI.BOOL], \
recvbuf=[Segment_Mean_Median[gtx].mask, MPI.BOOL], op=MPI.LAND)
Distributed_Mean_Median = None
# along-track X coordinates of segment fit
Segment_X_atc[gtx] = np.ma.zeros((n_seg),fill_value=fill_value)
Segment_X_atc[gtx].mask = np.ones((n_seg),dtype=bool)
comm.Allreduce(sendbuf=[Distributed_X_atc.data, MPI.DOUBLE], \
recvbuf=[Segment_X_atc[gtx].data, MPI.DOUBLE], op=MPI.SUM)
comm.Allreduce(sendbuf=[Distributed_X_atc.mask, MPI.BOOL], \
recvbuf=[Segment_X_atc[gtx].mask, MPI.BOOL], op=MPI.LAND)
Distributed_X_atc = None
# along-track X coordinate spread of points used in segment fit
Segment_X_spread[gtx] = np.ma.zeros((n_seg),fill_value=fill_value)
Segment_X_spread[gtx].mask = np.ones((n_seg),dtype=bool)
comm.Allreduce(sendbuf=[Distributed_X_spread.data, MPI.DOUBLE], \
recvbuf=[Segment_X_spread[gtx].data, MPI.DOUBLE], op=MPI.SUM)
comm.Allreduce(sendbuf=[Distributed_X_spread.mask, MPI.BOOL], \
recvbuf=[Segment_X_spread[gtx].mask, MPI.BOOL], op=MPI.LAND)
Distributed_X_spread = None
# along-track Y coordinates of segment fit
Segment_Y_atc[gtx] = np.ma.zeros((n_seg),fill_value=fill_value)
Segment_Y_atc[gtx].mask = np.ones((n_seg),dtype=bool)
comm.Allreduce(sendbuf=[Distributed_Y_atc.data, MPI.DOUBLE], \
recvbuf=[Segment_Y_atc[gtx].data, MPI.DOUBLE], op=MPI.SUM)
comm.Allreduce(sendbuf=[Distributed_Y_atc.mask, MPI.BOOL], \
recvbuf=[Segment_Y_atc[gtx].mask, MPI.BOOL], op=MPI.LAND)
Distributed_Y_atc = None
# longitude of fit photons
Segment_Longitude[gtx] = np.zeros((n_seg),fill_value=fill_value)
comm.Allreduce(sendbuf=[Distributed_Longitude, MPI.DOUBLE], \
recvbuf=[Segment_Longitude[gtx], MPI.DOUBLE], op=MPI.SUM)
# latitude of fit photons