-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathtag.py
238 lines (183 loc) · 8.81 KB
/
tag.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
"""
tag.py - TOMS Add Glint script
Add glint angle to TOMS (Total Ozone Mapping Spectrometer) HDF5 files
Usage:
find TOMSEPL2/1996/07/ -name '*.he5' | grep -v '_glint\.he5$' | sort | xargs python tag.py
Dependencies:
h5py; numpy; PySPICE; PDSImage; listfauxfile
- the latter two are provided in this repository
"""
import sys
import h5py
import spice
import numpy
import PDSImage
import listfauxfile
if __name__=="__main__" and sys.argv[1:]:
### Conversion from radans to degrees and back
rpd,dpr = spice.rpd(),spice.dpr()
### Load this python script as a SPICE kernel
spice.furnsh(__file__)
"""
Meta-kernel:
The files lister here must exist at the location specified
\begindata
KERNELS_TO_LOAD = (
'naif0011.tls'
'pck00010.tpc'
'de432s.bsp'
)
\begintext
"""
### Earth equatorial and polar radii; calculate flattening
re = spice.gdpool('BODY399_RADII',0,1)[1]
rp = spice.gdpool('BODY399_RADII',2,1)[1]
flattening = (re-rp) / re
######################################################################
def glintangle(scAltKm,scLonDeg,scLatDeg,surfLonDeg,surfLatDeg,uvEarth2Sun):
"""Function to calculate glint angle
Inputs:
S/C altitude, Longitude, Latitude
Surface point Longitude, Latitude (altitude is zero by definition)
Ephemeris time, s past J2000 EPOCH
"""
toms_alt_km = scAltKm
toms_lon_rad = rpd * scLonDeg ### Scale degrees to radians
toms_lat_rad = rpd * scLatDeg ### Scale degrees to radians
surf_lon_rad = rpd * surfLonDeg
surf_lat_rad = rpd * surfLatDeg
### Convert TOMS and surface points geodetic positions to ECEF vectors
toms_vec = spice.georec(toms_lon_rad,toms_lat_rad,toms_alt_km,re,flattening)
surfPoint = spice.georec(surf_lon_rad,surf_lat_rad,0.,re,flattening)
### Get normal at surface point, ECEF unit vector
uvSurfNormal = spice.surfnm(re,re,rp,surfPoint)
### Get cosine of incidence angle (angle between Sun & surface normal vector)
mu = spice.vdot(uvEarth2Sun,uvSurfNormal)
if mu <= 0.: return 999.9
uvReflect = spice.vsub(spice.vscl(2.*mu,uvSurfNormal),uvEarth2Sun) ### Get specular reflection vector
### Return glint angle
return dpr * spice.vsep(uvReflect,spice.vsub(toms_vec,surfPoint))
### End of glint angle calculation
######################################################################
######################################################################
for fnInput in (__name__=="__main__" and sys.argv[1:] or []):
if fnInput[-4:] != '.he5':
sys.stderr.write("### Skipping file '%s'; does not end in '.he5'\n" % (fnInput,))
continue
if fnInput[-10:] == '_glint.he5':
sys.stderr.write("### Skipping file '%s'; ends in '_glint.he5'\n" % (fnInput,))
continue
### Build output filename
fnToks = fnInput.split('.')
fnToks[-2] += '_glint'
fnOutput = '.'.join(fnToks)
try:
fOutput = h5py.File(fnOutput,'w-')
except:
sys.stderr.write("### Skipping creation of glint file '%s'; it already exists or is otherwise un-writeable\n" % (fnOutput,))
continue
### Read HDF5 file
fInput=h5py.File(fnInput,'r')
### Get the inputs from the input HDF5 file
getThem = lambda s: numpy.array(fInput[s])
tomsOffsetTimes,surfLatDegs,surfLonDegs,scAltKms,scLatDegs,scLonDegs = map(getThem
, [s.strip() for s in """
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/Time
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/Latitude
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/Longitude
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/SpacecraftAltitude
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/SpacecraftLatitude
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/SpacecraftLongitude
""".strip().split('\n')]
)
### Convert .../CoreMetadata to string, parse as ODL using PDSImage.py
faux=listfauxfile.FILE(numpy.array(fInput['HDFEOS INFORMATION/CoreMetadata']).tostring().replace('\0',''))
pdsi=PDSImage.PDSImage()
### dcmd ocontains recursive dictionaries
dcmd=pdsi._downone(faux)
### Extract start time from dcmd (CoreMetadata) as ISO UTC string
firstTomsUTC = '%sT%s' % (dcmd['INVENTORYMETADATA']['RANGEDATETIME']['RANGEBEGINNINGDATE']['VALUE']
,dcmd['INVENTORYMETADATA']['RANGEDATETIME']['RANGEBEGINNINGTIME']['VALUE']
,)
### Combine start time and /Time offset to get the TOMS data zero epoch as seconds past J2000 epoch
### N.B. is typically (always?) equivalent to 1993-01-01T00:00:00 UTC
tomsZeroEpoch = spice.utc2et(firstTomsUTC) - tomsOffsetTimes[0]
print('%s(Epoch) %s(Start) %s' % (spice.et2utc(tomsZeroEpoch,'ISOC',1),firstTomsUTC,fnInput,))
### Iterators
twoDshape = surfLonDegs.shape
rowIter,pixelIter = map(xrange,twoDshape)
glintAngles = numpy.zeros(twoDshape,dtype=numpy.float32)
for row in rowIter:
### Earth-Sun vector, convert to unit vector
earth2Sun = spice.spkezr('SUN',tomsZeroEpoch+tomsOffsetTimes[row],'IAU_EARTH','LT+S','EARTH')[0][:3]
uvEarth2Sun = spice.vhat(earth2Sun)
scAltKm = scAltKms[row]
scLonDeg = scLonDegs[row]
scLatDeg = scLatDegs[row]
for pixel in pixelIter:
glintAngles[row,pixel] = glintangle(scAltKm
,scLonDeg
,scLatDeg
,surfLonDegs[row,pixel]
,surfLatDegs[row,pixel]
,uvEarth2Sun
)
### Get top level groups for copying
topGroups = []
def addtopgroup(name):
if name.find('/')==-1: topGroups.append(name)
fInput.visit(addtopgroup)
for topGroup in topGroups: fOutput.copy( fInput[topGroup], topGroup)
fInput.close()
fOutput.create_dataset('HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/GlintAngle',data=glintAngles)
fOutput.close()
"""
HDFEOS
HDFEOS/ADDITIONAL
HDFEOS/ADDITIONAL/FILE_ATTRIBUTES
HDFEOS/SWATHS
HDFEOS/SWATHS/EP TOMS Column Amount O3
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/APrioriLayerO3
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/AlgorithmFlags
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/CloudFraction
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/CloudTopPressure
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/ColumnAmountO3
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/LayerEfficiency
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/MeasurementQualityFlags
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/NValue
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/O3BelowCloud
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/QualityFlags
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/Reflectivity331
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/Reflectivity360
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/Residual
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/ResidualStep1
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/ResidualStep2
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/SO2index
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/Sensitivity
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/StepOneO3
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/StepTwoO3
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/TerrainPressure
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/UVAerosolIndex
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/Wavelength
HDFEOS/SWATHS/EP TOMS Column Amount O3/Data Fields/dN_dR
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/GroundPixelQualityFlags
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/Latitude
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/Longitude
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/RelativeAzimuthAngle
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/SecondsInDay
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/SolarAzimuthAngle
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/SolarZenithAngle
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/SpacecraftAltitude
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/SpacecraftLatitude
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/SpacecraftLongitude
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/TerrainHeight
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/Time
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/ViewingAzimuthAngle
HDFEOS/SWATHS/EP TOMS Column Amount O3/Geolocation Fields/ViewingZenithAngle
HDFEOS INFORMATION
HDFEOS INFORMATION/ArchivedMetadata
HDFEOS INFORMATION/CoreMetadata
HDFEOS INFORMATION/StructMetadata.0
"""