1
|
|
|
#!/usr/bin/env python |
2
|
|
|
# vim: set fileencoding=utf-8 |
3
|
|
|
"""SCIAMACHY level 2 data post processing |
4
|
|
|
|
5
|
|
|
Main script for SCIAMACHY orbital retrieval post processing |
6
|
|
|
and data combining (to netcdf). |
7
|
|
|
|
8
|
|
|
Copyright (c) 2018 Stefan Bender |
9
|
|
|
|
10
|
|
|
This file is part of sciapy. |
11
|
|
|
sciapy is free software: you can redistribute it or modify it |
12
|
|
|
under the terms of the GNU General Public License as published by |
13
|
|
|
the Free Software Foundation, version 2. |
14
|
|
|
See accompanying LICENSE file or http://www.gnu.org/licenses/gpl-2.0.html. |
15
|
|
|
""" |
16
|
|
|
|
17
|
|
|
from __future__ import absolute_import, division, print_function |
18
|
|
|
|
19
|
|
|
import glob |
20
|
|
|
import os |
21
|
|
|
import argparse as ap |
22
|
|
|
import datetime as dt |
23
|
|
|
import logging |
24
|
|
|
from pkg_resources import resource_filename |
25
|
|
|
|
26
|
|
|
import numpy as np |
27
|
|
|
import pandas as pd |
28
|
|
|
import xarray as xr |
29
|
|
|
from scipy.interpolate import interp1d |
30
|
|
|
|
31
|
|
|
from astropy import units |
32
|
|
|
from astropy.time import Time |
33
|
|
|
import astropy.coordinates as coord |
34
|
|
|
|
35
|
|
|
import sciapy.level1c as sl |
36
|
|
|
from sciapy.level2 import density_pp as sd |
37
|
|
|
from sciapy.level2 import scia_akm as sa |
38
|
|
|
from sciapy.level2.igrf import gmag_igrf |
39
|
|
|
from sciapy.level2.aacgm2005 import gmag_aacgm2005 |
40
|
|
|
try: |
41
|
|
|
from nrlmsise00 import msise_flat as msise |
42
|
|
|
except ImportError: |
43
|
|
|
msise = None |
44
|
|
|
try: |
45
|
|
|
from sciapy.level2.noem import noem_cpp |
46
|
|
|
except ImportError: |
47
|
|
|
noem_cpp = None |
48
|
|
|
|
49
|
|
|
F107_FILE = resource_filename("sciapy", "data/indices/f107_noontime_flux_obs.txt") |
50
|
|
|
F107A_FILE = resource_filename("sciapy", "data/indices/f107a_noontime_flux_obs.txt") |
51
|
|
|
AP_FILE = resource_filename("sciapy", "data/indices/spidr_ap_2000-2012.dat") |
52
|
|
|
F107_ADJ_FILE = resource_filename("sciapy", "data/indices/spidr_f107_2000-2012.dat") |
53
|
|
|
KP_FILE = resource_filename("sciapy", "data/indices/spidr_kp_2000-2012.dat") |
54
|
|
|
|
55
|
|
|
PHI_FAC = 11.91 |
56
|
|
|
LST_FAC = -0.62 |
57
|
|
|
|
58
|
|
|
|
59
|
|
View Code Duplication |
def solar_zenith_angle(alt, lat, lon, time): |
|
|
|
|
60
|
|
|
atime = Time(time) |
61
|
|
|
loc = coord.EarthLocation.from_geodetic( |
62
|
|
|
height=alt * units.km, |
63
|
|
|
lat=lat * units.deg, |
64
|
|
|
lon=lon * units.deg, |
65
|
|
|
) |
66
|
|
|
altaz = coord.AltAz(location=loc, obstime=atime) |
67
|
|
|
sun = coord.get_sun(atime) |
68
|
|
|
return sun.transform_to(altaz).zen.value |
69
|
|
|
|
70
|
|
|
|
71
|
|
View Code Duplication |
def read_spectra(year, orbit, spec_base=None, skip_upleg=True): |
|
|
|
|
72
|
|
|
"""Read and examine SCIAMACHY orbit spectra |
73
|
|
|
|
74
|
|
|
Reads the limb spactra and extracts the dates, times, latitudes, |
75
|
|
|
longitudes to be used to re-assess the retrieved geolocations. |
76
|
|
|
|
77
|
|
|
Parameters |
78
|
|
|
---------- |
79
|
|
|
year: int |
80
|
|
|
The measurement year to select the corresponding subdir |
81
|
|
|
below `spec_base` (see below). |
82
|
|
|
orbit: int |
83
|
|
|
SCIAMACHY/Envisat orbit number of the spectra. |
84
|
|
|
spec_base: str, optional |
85
|
|
|
The root path to the level 1c spectra. Uses the current |
86
|
|
|
dir if not set or set to `None` (default). |
87
|
|
|
skip_upleg: bool, optional |
88
|
|
|
Skip upleg limb scans, i.e. night time scans. For NO retrievals, |
89
|
|
|
those are not used and should be not used here. |
90
|
|
|
Default: True |
91
|
|
|
|
92
|
|
|
Returns |
93
|
|
|
------- |
94
|
|
|
(dts, times, lats, lons, mlsts, alsts, eotcorr) |
95
|
|
|
""" |
96
|
|
|
fail = (None,) * 7 |
97
|
|
|
if spec_base is None: |
98
|
|
|
spec_base = os.curdir |
99
|
|
|
spec_path = os.path.join(spec_base, "{0}".format(year)) |
100
|
|
|
spec_path2 = os.path.join(spec_base, "{0}".format(int(year) + 1)) |
101
|
|
|
logging.debug("spec_path: %s", spec_path) |
102
|
|
|
logging.debug("spec_path2: %s", spec_path) |
103
|
|
|
if not (os.path.isdir(spec_path) or os.path.isdir(spec_path2)): |
104
|
|
|
return fail |
105
|
|
|
|
106
|
|
|
# the star stands for the (optional) date subdir |
107
|
|
|
# to find all spectra for the orbit |
108
|
|
|
spfiles = glob.glob( |
109
|
|
|
'{0}/*/SCIA_limb_*_1_0_{1:05d}.dat.l_mpl_binary' |
110
|
|
|
.format(spec_path, orbit)) |
111
|
|
|
# sometimes for whatever reason the orbit ends up in the wrong year subdir |
112
|
|
|
# looks in the subdir for the following year as well. |
113
|
|
|
spfiles += glob.glob( |
114
|
|
|
'{0}/*/SCIA_limb_*_1_0_{1:05d}.dat.l_mpl_binary' |
115
|
|
|
.format(spec_path2, orbit)) |
116
|
|
|
if len(spfiles) < 2: |
117
|
|
|
return fail |
118
|
|
|
|
119
|
|
|
dts = [] |
120
|
|
|
times = [] |
121
|
|
|
lats = [] |
122
|
|
|
lons = [] |
123
|
|
|
mlsts = [] |
124
|
|
|
alsts = [] |
125
|
|
|
|
126
|
|
|
sls = sl.scia_limb_scan() |
127
|
|
|
for f in sorted(spfiles): |
128
|
|
|
sls.read_from_file(f) |
129
|
|
|
# copy the values from the l1c file |
130
|
|
|
lat, lon = sls.cent_lat_lon[:2] |
131
|
|
|
mlst, alst, eotcorr = sls.local_solar_time(False) |
132
|
|
|
tp_lats = sls.limb_data.tp_lat |
133
|
|
|
date = sls.date |
134
|
|
|
# debug output if requested |
135
|
|
|
logging.debug("file: %s", f) |
136
|
|
|
logging.debug("lat: %s, lon: %s", lat, lon) |
137
|
|
|
logging.debug("mlst: %s, alst: %s, eotcorr: %s", mlst, alst, eotcorr) |
138
|
|
|
logging.debug("tp_lats: %s", tp_lats) |
139
|
|
|
logging.debug("date: %s", date) |
140
|
|
|
if skip_upleg and ((tp_lats[1] - tp_lats[-2]) < 0.5): |
141
|
|
|
# Exclude non-downleg measurements where the latitude |
142
|
|
|
# of the last real tangent point (the last is dark sky) |
143
|
|
|
# is larger than or too close to the first latitude. |
144
|
|
|
# Requires an (empirical) separation of +0.5 degree. |
145
|
|
|
logging.debug("excluding upleg point at: %s, %s", lat, lon) |
146
|
|
|
continue |
147
|
|
|
dtdate = pd.to_datetime(dt.datetime(*date), utc=True) |
148
|
|
|
time_hour = dtdate.hour + dtdate.minute / 60.0 + dtdate.second / 3600.0 |
149
|
|
|
logging.debug("mean lst: %s, apparent lst: %s, EoT: %s", mlst, alst, eotcorr) |
150
|
|
|
dts.append(dtdate) |
151
|
|
|
times.append(time_hour) |
152
|
|
|
lats.append(lat) |
153
|
|
|
lons.append(lon) |
154
|
|
|
mlsts.append(mlst) |
155
|
|
|
alsts.append(alst) |
156
|
|
|
|
157
|
|
|
if len(lats) < 2: |
158
|
|
|
# interpolation will fail with less than 2 points |
159
|
|
|
return fail |
160
|
|
|
|
161
|
|
|
return (np.asarray(dts), |
162
|
|
|
np.asarray(times), |
163
|
|
|
np.asarray(lats), |
164
|
|
|
np.asarray(lons), |
165
|
|
|
np.asarray(mlsts), |
166
|
|
|
np.asarray(alsts), eotcorr) |
|
|
|
|
167
|
|
|
|
168
|
|
|
|
169
|
|
|
class _circ_interp(object): |
170
|
|
|
"""Interpolation on a circle""" |
171
|
|
|
def __init__(self, x, y, **kw): |
172
|
|
|
self.c_intpf = interp1d(x, np.cos(y), **kw) |
173
|
|
|
self.s_intpf = interp1d(x, np.sin(y), **kw) |
174
|
|
|
|
175
|
|
|
def __call__(self, x): |
176
|
|
|
return np.arctan2(self.s_intpf(x), self.c_intpf(x)) |
177
|
|
|
|
178
|
|
|
|
179
|
|
View Code Duplication |
def process_orbit( |
|
|
|
|
180
|
|
|
orbit, |
181
|
|
|
ref_date="1950-01-01", |
182
|
|
|
dens_path=None, |
183
|
|
|
spec_base=None, |
184
|
|
|
use_msis=True, |
185
|
|
|
): |
186
|
|
|
"""Post process retrieved SCIAMACHY orbit |
187
|
|
|
|
188
|
|
|
Parameters |
189
|
|
|
---------- |
190
|
|
|
orbit: int |
191
|
|
|
SCIAMACHY/Envisat orbit number of the results to process. |
192
|
|
|
ref_date: str, optional |
193
|
|
|
Base date to calculate the relative days from, |
194
|
|
|
of the format "%Y-%m-%d". Default: 1950-01-01 |
195
|
|
|
dens_path: str, optional |
196
|
|
|
The path to the level 2 data. Uses the current |
197
|
|
|
dir if not set or set to `None` (default). |
198
|
|
|
spec_base: str, optional |
199
|
|
|
The root path to the level 1c spectra. Uses the current |
200
|
|
|
dir if not set or set to `None` (default). |
201
|
|
|
|
202
|
|
|
Returns |
203
|
|
|
------- |
204
|
|
|
(dts0, time0, lst0, lon0, sdd): tuple |
205
|
|
|
dts0 - days since ref_date at equator crossing (float) |
206
|
|
|
time0 - utc hour into the day at equator crossing (float) |
207
|
|
|
lst0 - apparent local solar time at the equator (float) |
208
|
|
|
lon0 - longitude of the equator crossing (float) |
209
|
|
|
sdd - `scia_density_pp` instance of the post-processed data |
210
|
|
|
""" |
211
|
|
|
def _read_gm(fname): |
212
|
|
|
return dict(np.genfromtxt(fname, usecols=[0, 2], dtype=None)) |
213
|
|
|
|
214
|
|
|
fail = (None,) * 5 |
215
|
|
|
logging.debug("processing orbit: %s", orbit) |
216
|
|
|
dtrefdate = pd.to_datetime(ref_date, format="%Y-%m-%d", utc=True) |
217
|
|
|
|
218
|
|
|
dfiles = glob.glob( |
219
|
|
|
"{0}/000NO_orbit_{1:05d}_*_Dichten.txt" |
220
|
|
|
.format(dens_path, orbit)) |
221
|
|
|
if len(dfiles) < 1: |
222
|
|
|
return fail |
223
|
|
|
logging.debug("dfiles: %s", dfiles) |
224
|
|
|
logging.debug("splits: %s", [fn.split('/') for fn in dfiles]) |
225
|
|
|
ddict = dict([(fn, (fn.split('/')[-3:-1] + fn.split('/')[-1].split('_')[3:4])) |
226
|
|
|
for fn in dfiles]) |
227
|
|
|
logging.debug("ddict: %s", ddict) |
228
|
|
|
year = ddict[sorted(ddict.keys())[0]][-1][:4] |
229
|
|
|
logging.debug("year: %s", year) |
230
|
|
|
|
231
|
|
|
dts, times, lats, lons, mlsts, alsts, eotcorr = \ |
232
|
|
|
read_spectra(year, orbit, spec_base) |
233
|
|
|
if dts is None: |
234
|
|
|
# return early if reading the spectra failed |
235
|
|
|
return fail |
236
|
|
|
|
237
|
|
|
dts = pd.to_datetime(dts) - dtrefdate |
238
|
|
|
dts = np.array([dtd.days + dtd.seconds / 86400. for dtd in dts]) |
239
|
|
|
logging.debug("lats: %s, lons: %s, times: %s", lats, lons, times) |
240
|
|
|
|
241
|
|
|
sdd = sd.scia_densities_pp(ref_date=ref_date) |
242
|
|
|
sdd.read_from_file(dfiles[0]) |
243
|
|
|
logging.debug("density lats: %s, lons: %s", sdd.lats, sdd.lons) |
244
|
|
|
|
245
|
|
|
# Re-interpolates the location (longitude) and times from the |
246
|
|
|
# limb scan spectra files along the orbit to determine the values |
247
|
|
|
# at the Equator and to fill in possibly missing data. |
248
|
|
|
# |
249
|
|
|
# y values are unit circle angles in radians (0 < φ < 2π or -π < φ < π) |
250
|
|
|
# longitudes |
251
|
|
|
lons_intpf = _circ_interp( |
252
|
|
|
lats[::-1], np.radians(lons[::-1]), |
253
|
|
|
fill_value="extrapolate", |
254
|
|
|
) |
255
|
|
|
# apparent local solar time (EoT corrected) |
256
|
|
|
lst_intpf = _circ_interp( |
257
|
|
|
lats[::-1], np.pi / 12. * alsts[::-1], |
258
|
|
|
fill_value="extrapolate", |
259
|
|
|
) |
260
|
|
|
# mean local solar time |
261
|
|
|
mst_intpf = _circ_interp( |
262
|
|
|
lats[::-1], np.pi / 12. * mlsts[::-1], |
263
|
|
|
fill_value="extrapolate", |
264
|
|
|
) |
265
|
|
|
# utc time (day) |
266
|
|
|
time_intpf = _circ_interp( |
267
|
|
|
lats[::-1], np.pi / 12. * times[::-1], |
268
|
|
|
fill_value="extrapolate", |
269
|
|
|
) |
270
|
|
|
# datetime |
271
|
|
|
dts_retr_interpf = interp1d(lats[::-1], dts[::-1], fill_value="extrapolate") |
272
|
|
|
|
273
|
|
|
# equator values |
274
|
|
|
lon0 = np.degrees(lons_intpf(0.)) % 360. |
275
|
|
|
lst0 = (lst_intpf(0.) * 12. / np.pi) % 24. |
276
|
|
|
mst0 = (mst_intpf(0.) * 12. / np.pi) % 24. |
277
|
|
|
time0 = (time_intpf(0.) * 12. / np.pi) % 24. |
278
|
|
|
dts_retr_interp0 = dts_retr_interpf(0.) |
279
|
|
|
logging.debug("utc day at equator: %s", dts_retr_interp0) |
280
|
|
|
logging.debug("mean LST at equator: %s, apparent LST at equator: %s", mst0, lst0) |
281
|
|
|
|
282
|
|
|
sdd.utchour = (time_intpf(sdd.lats) * 12. / np.pi) % 24. |
283
|
|
|
sdd.utcdays = dts_retr_interpf(sdd.lats) |
284
|
|
|
|
285
|
|
|
if sdd.lons is None: |
286
|
|
|
# recalculate the longitudes |
287
|
|
|
# estimate the equatorial longitude from the |
288
|
|
|
# limb scan latitudes and longitudes |
289
|
|
|
lon0s_tp = lons - PHI_FAC * np.tan(np.radians(lats)) |
290
|
|
|
clon0s_tp = np.cos(np.radians(lon0s_tp)) |
291
|
|
|
slon0s_tp = np.sin(np.radians(lon0s_tp)) |
292
|
|
|
lon0_tp = np.arctan2(np.sum(slon0s_tp[1:-1]), np.sum(clon0s_tp[1:-1])) |
293
|
|
|
lon0_tp = np.degrees((lon0_tp + 2. * np.pi) % (2. * np.pi)) |
294
|
|
|
logging.info("lon0: %s", lon0) |
295
|
|
|
logging.info("lon0 tp: %s", lon0_tp) |
296
|
|
|
# interpolate to the retrieval latitudes |
297
|
|
|
tg_retr_lats = np.tan(np.radians(sdd.lats)) |
298
|
|
|
calc_lons = (tg_retr_lats * PHI_FAC + lon0) % 360. |
299
|
|
|
calc_lons_tp = (tg_retr_lats * PHI_FAC + lon0_tp) % 360. |
300
|
|
|
sdd.lons = calc_lons_tp |
301
|
|
|
logging.debug("(calculated) retrieval lons: %s, %s", |
302
|
|
|
calc_lons, calc_lons_tp) |
303
|
|
|
else: |
304
|
|
|
# sdd.lons = sdd.lons % 360. |
305
|
|
|
logging.debug("(original) retrieval lons: %s", sdd.lons) |
306
|
|
|
|
307
|
|
|
sdd.mst = (sdd.utchour + sdd.lons / 15.) % 24. |
308
|
|
|
sdd.lst = sdd.mst + eotcorr / 60. |
309
|
|
|
mean_alt_km = sdd.alts.mean() |
310
|
|
|
mean_alt_m = mean_alt_km * 1000. |
311
|
|
|
|
312
|
|
|
dt_date_this = dt.timedelta(np.asscalar(dts_retr_interp0)) + dtrefdate |
313
|
|
|
logging.info("date: %s", dt_date_this) |
314
|
|
|
# caclulate geomagnetic coordinates |
315
|
|
|
sdd.gmlats, sdd.gmlons = gmag_igrf(dt_date_this, sdd.lats, sdd.lons, alt=mean_alt_km) |
316
|
|
|
logging.debug("geomag. lats: %s, lons: %s", sdd.gmlats, sdd.gmlons) |
317
|
|
|
sdd.aacgmgmlats, sdd.aacgmgmlons = gmag_aacgm2005(sdd.lats, sdd.lons) |
318
|
|
|
logging.debug("aacgm geomag. lats: %s, lons: %s", |
319
|
|
|
sdd.aacgmgmlats, sdd.aacgmgmlons) |
320
|
|
|
|
321
|
|
|
# current day for MSIS input |
322
|
|
|
f107_data = _read_gm(F107_FILE) |
323
|
|
|
f107a_data = _read_gm(F107A_FILE) |
324
|
|
|
ap_data = _read_gm(AP_FILE) |
325
|
|
|
msis_dtdate = dt.timedelta(np.asscalar(dts_retr_interp0)) + dtrefdate |
326
|
|
|
msis_dtdate1 = msis_dtdate - dt.timedelta(days=1) |
327
|
|
|
msis_date = msis_dtdate.strftime("%Y-%m-%d").encode() |
328
|
|
|
msis_date1 = msis_dtdate1.strftime("%Y-%m-%d").encode() |
329
|
|
|
msis_f107 = f107_data[msis_date1] |
330
|
|
|
msis_f107a = f107a_data[msis_date] |
331
|
|
|
msis_ap = ap_data[msis_date] |
332
|
|
|
logging.debug("MSIS date: %s, f10.7a: %s, f10.7: %s, ap: %s", |
333
|
|
|
msis_date, msis_f107a, msis_f107, msis_ap) |
334
|
|
|
|
335
|
|
|
# previous day for NOEM input |
336
|
|
|
f107_adj = _read_gm(F107_ADJ_FILE) |
337
|
|
|
kp_data = _read_gm(KP_FILE) |
338
|
|
|
noem_dtdate = dt.timedelta(np.asscalar(dts_retr_interp0) - 1) + dtrefdate |
339
|
|
|
noem_date = noem_dtdate.strftime("%Y-%m-%d").encode() |
340
|
|
|
noem_f107 = f107_adj[noem_date] |
341
|
|
|
noem_kp = kp_data[noem_date] |
342
|
|
|
logging.debug("NOEM date: %s, f10.7: %s, kp: %s", |
343
|
|
|
noem_date, noem_f107, noem_kp) |
344
|
|
|
|
345
|
|
|
if sdd.noem_no is None: |
346
|
|
|
sdd.noem_no = np.zeros_like(sdd.densities) |
347
|
|
|
if sdd.temperature is None and msise is None: |
348
|
|
|
sdd.temperature = np.full_like(sdd.densities, np.nan) |
349
|
|
|
if sdd.sza is None: |
350
|
|
|
sdd.sza = np.zeros_like(sdd.lats) |
351
|
|
|
if sdd.akdiag is None: |
352
|
|
|
sdd.akdiag = np.zeros_like(sdd.densities) |
353
|
|
|
akm_filename = glob.glob( |
354
|
|
|
"{0}/000NO_orbit_{1:05d}_*_AKM*" |
355
|
|
|
.format(dens_path, orbit))[0] |
356
|
|
|
logging.debug("ak file: %s", akm_filename) |
357
|
|
|
ak = sa.read_akm(akm_filename, sdd.nalt, sdd.nlat) |
358
|
|
|
logging.debug("ak data: %s", ak) |
359
|
|
|
sdd.akdiag = ak.diagonal(axis1=1, axis2=3).diagonal(axis1=0, axis2=1) |
360
|
|
|
|
361
|
|
|
if msise is not None: |
362
|
|
|
if sdd.temperature is None or use_msis: |
363
|
|
|
_msis_d_t = msise( |
364
|
|
|
msis_dtdate, |
365
|
|
|
sdd.alts[None, :], sdd.lats[:, None], sdd.lons[:, None] % 360., |
366
|
|
|
msis_f107a, msis_f107, msis_ap, |
367
|
|
|
lst=sdd.lst[:, None], |
368
|
|
|
) |
369
|
|
|
sdd.temperature = _msis_d_t[:, :, -1] |
370
|
|
|
if use_msis: |
371
|
|
|
sdd.dens_tot = np.sum(_msis_d_t[:, :, np.r_[:5, 6:9]], axis=2) |
|
|
|
|
372
|
|
|
for i, lat in enumerate(sdd.lats): |
373
|
|
|
if noem_cpp is not None: |
374
|
|
|
sdd.noem_no[i] = noem_cpp(noem_date.decode(), sdd.alts, |
375
|
|
|
[lat], [sdd.lons[i]], noem_f107, noem_kp)[:] |
376
|
|
|
else: |
377
|
|
|
sdd.noem_no[i][:] = np.nan |
378
|
|
|
sdd.sza[:] = solar_zenith_angle( |
379
|
|
|
mean_alt_km, |
380
|
|
|
sdd.lats, sdd.lons, |
381
|
|
|
pd.to_timedelta(sdd.utcdays, unit="days") + dtrefdate, |
382
|
|
|
) |
383
|
|
|
sdd.vmr = sdd.densities / sdd.dens_tot * 1.e9 # ppb |
384
|
|
|
return dts_retr_interp0, time0, lst0, lon0, sdd |
385
|
|
|
|
386
|
|
|
|
387
|
|
View Code Duplication |
def get_orbits_from_date(date, mlt=False, path=None, L2_version="v6.2"): |
|
|
|
|
388
|
|
|
"""Find SCIAMACHY orbits with retrieved data at a date |
389
|
|
|
|
390
|
|
|
Parameters |
391
|
|
|
---------- |
392
|
|
|
date: str |
393
|
|
|
The date in the format "%Y-%m-%d". |
394
|
|
|
mlt: bool, optional |
395
|
|
|
Look for MLT mode data instead of nominal mode data. |
396
|
|
|
Increases the heuristics to find all MLT orbits. |
397
|
|
|
path: str, optional |
398
|
|
|
The path to the level 2 data. If `None` tries to infer |
399
|
|
|
the data directory from the L2 version using |
400
|
|
|
'./*<L2_version>'. Default: None |
401
|
|
|
|
402
|
|
|
Returns |
403
|
|
|
------- |
404
|
|
|
orbits: list |
405
|
|
|
List of found orbits with retrieved data files |
406
|
|
|
""" |
407
|
|
|
logging.debug("pre-processing: %s", date) |
408
|
|
|
if path is None: |
409
|
|
|
density_base = os.curdir |
410
|
|
|
path = "{0}/*{1}".format(density_base, L2_version) |
411
|
|
|
logging.debug("path: %s", path) |
412
|
|
|
|
413
|
|
|
dfiles = glob.glob("{0}/000NO_orbit_*_{1}_Dichten.txt".format( |
414
|
|
|
path, date.replace("-", ""))) |
415
|
|
|
orbits = sorted([int(os.path.basename(df).split('_')[2]) for df in dfiles]) |
416
|
|
|
if mlt: |
417
|
|
|
orbits.append(orbits[-1] + 1) |
418
|
|
|
return orbits |
419
|
|
|
|
420
|
|
|
|
421
|
|
View Code Duplication |
def combine_orbit_data(orbits, |
|
|
|
|
422
|
|
|
ref_date="1950-01-01", |
423
|
|
|
L2_version="v6.2", file_version="2.3", |
424
|
|
|
dens_path=None, spec_base=None, |
425
|
|
|
use_xarray=False, save_nc=False): |
426
|
|
|
"""Combine post-processed SCIAMACHY retrieved orbit data |
427
|
|
|
|
428
|
|
|
Parameters |
429
|
|
|
---------- |
430
|
|
|
orbits: list |
431
|
|
|
List of SCIAMACHY/Envisat orbit numbers to process. |
432
|
|
|
ref_date: str, optional |
433
|
|
|
Base date to calculate the relative days from, |
434
|
|
|
of the format "%Y-%m-%d". Default: 1950-01-01 |
435
|
|
|
L2_version: str, optional |
436
|
|
|
SCIAMACHY level 2 data version to process |
437
|
|
|
file_version: str, optional |
438
|
|
|
Postprocessing format version of the output data |
439
|
|
|
dens_path: str, optional |
440
|
|
|
The path to the level 2 data. If `None` tries to infer |
441
|
|
|
the data directory from the L2 version looking for anything |
442
|
|
|
in the current directory that ends in <L2_version>: './*<L2_version>'. |
443
|
|
|
Default: None |
444
|
|
|
spec_base: str, optional |
445
|
|
|
The root path to the level 1c spectra. Uses the current |
446
|
|
|
dir if not set or set to `None` (default). |
447
|
|
|
use_xarray: bool, optional |
448
|
|
|
Uses xarray (if available) to combine the orbital data. |
449
|
|
|
save_nc: bool, optional |
450
|
|
|
Save the intermediate orbit data sets to netcdf files |
451
|
|
|
for debugging. |
452
|
|
|
|
453
|
|
|
Returns |
454
|
|
|
------- |
455
|
|
|
(sdday, sdday_ds): tuple |
456
|
|
|
`sdday` contains the combined data as a `scia_density_day` instance, |
457
|
|
|
`sdday_ds` contains the same data as a `xarray.Dataset`. |
458
|
|
|
""" |
459
|
|
|
if dens_path is None: |
460
|
|
|
# try some heuristics |
461
|
|
|
density_base = os.curdir |
462
|
|
|
dens_path = "{0}/*{1}".format(density_base, L2_version) |
463
|
|
|
|
464
|
|
|
sdday = sd.scia_density_day(ref_date=ref_date) |
465
|
|
|
sddayl = [] |
466
|
|
|
sdday_ds = None |
467
|
|
|
for orbit in sorted(orbits): |
468
|
|
|
dateo, timeo, lsto, lono, sdens = process_orbit(orbit, |
469
|
|
|
ref_date=ref_date, dens_path=dens_path, spec_base=spec_base) |
470
|
|
|
logging.info( |
471
|
|
|
"orbit: %s, eq. date: %s, eq. hour: %s, eq. app. lst: %s, eq. lon: %s", |
472
|
|
|
orbit, dateo, timeo, lsto, lono |
473
|
|
|
) |
474
|
|
|
if sdens is not None: |
475
|
|
|
sdens.version = file_version |
476
|
|
|
sdens.data_version = L2_version |
477
|
|
|
sdday.append_data(dateo, orbit, timeo, sdens) |
478
|
|
|
if use_xarray: |
479
|
|
|
sd_xr = sdens.to_xarray(dateo, orbit) |
480
|
|
|
if sd_xr is not None: |
481
|
|
|
logging.debug("orbit %s dataset: %s", orbit, sd_xr) |
482
|
|
|
sddayl.append(sd_xr) |
483
|
|
|
if save_nc: |
484
|
|
|
sdens.write_to_netcdf(sdens.filename[:-3] + "nc") |
485
|
|
|
if use_xarray and sddayl: |
486
|
|
|
sdday_ds = xr.concat(sddayl, dim="time") |
487
|
|
|
return sdday, sdday_ds |
488
|
|
|
|
489
|
|
|
|
490
|
|
View Code Duplication |
def sddata_xr_set_attrs(sdday_xr, ref_date="1950-01-01", rename=True): |
|
|
|
|
491
|
|
|
"""Customize xarray Dataset variables and attributes |
492
|
|
|
|
493
|
|
|
Changes the variable names to match those exported from the |
494
|
|
|
`scia_density_day` class. |
495
|
|
|
|
496
|
|
|
Parameters |
497
|
|
|
---------- |
498
|
|
|
sdday_xr: xarray.Dataset instance |
499
|
|
|
ref_date: str, optional |
500
|
|
|
Base date to calculate the relative days from, |
501
|
|
|
of the format "%Y-%m-%d". Default: 1950-01-01 |
502
|
|
|
rename: bool, optional |
503
|
|
|
Rename the dataset variables to match the |
504
|
|
|
`scia_density_day` exported ones. |
505
|
|
|
""" |
506
|
|
|
if rename: |
507
|
|
|
sdday_xr = sdday_xr.rename({ |
508
|
|
|
# 2d vars |
509
|
|
|
"akm_diagonal": "NO_AKDIAG", |
510
|
|
|
"apriori": "NO_APRIORI", |
511
|
|
|
"density": "NO_DENS", |
512
|
|
|
"density_air": "MSIS_Dens", |
513
|
|
|
"error_meas": "NO_ERR", |
514
|
|
|
"error_tot": "NO_ETOT", |
515
|
|
|
"temperature": "MSIS_Temp", |
516
|
|
|
"NOEM_density": "NO_NOEM", |
517
|
|
|
"VMR": "NO_VMR", |
518
|
|
|
# 1d vars and dimensions |
519
|
|
|
"app_lst": "app_LST", |
520
|
|
|
"mean_lst": "mean_LST", |
521
|
|
|
"mean_sza": "mean_SZA", |
522
|
|
|
"utc_hour": "UTC", |
523
|
|
|
}) |
524
|
|
|
sdday_xr["NO_RSTD"] = 100 * np.abs(sdday_xr.NO_ERR / sdday_xr.NO_DENS) |
525
|
|
|
sdday_xr["NO_RSTD"].attrs = dict(units='%', |
526
|
|
|
long_name='NO relative standard deviation') |
527
|
|
|
# fix coordinate attributes |
528
|
|
|
sdday_xr["time"].attrs = dict(axis='T', standard_name='time', |
529
|
|
|
calendar='standard', long_name='equatorial crossing time', |
530
|
|
|
units="days since {0}".format( |
531
|
|
|
pd.to_datetime(ref_date).isoformat(sep=" "))) |
532
|
|
|
sdday_xr["orbit"].attrs = dict(axis='T', calendar='standard', |
533
|
|
|
long_name='SCIAMACHY/Envisat orbit number', units='1') |
534
|
|
|
sdday_xr["altitude"].attrs = dict(axis='Z', positive='up', |
535
|
|
|
long_name='altitude', standard_name='altitude', units='km') |
536
|
|
|
sdday_xr["latitude"].attrs = dict(axis='Y', long_name='latitude', |
537
|
|
|
standard_name='latitude', units='degrees_north') |
538
|
|
|
sdday_xr["longitude"].attrs = dict(long_name='longitude', |
539
|
|
|
standard_name='longitude', units='degrees_east') |
540
|
|
|
dateo = (pd.to_datetime( |
541
|
|
|
xr.conventions.decode_cf_variable("date", sdday_xr.time).data[0]) |
542
|
|
|
.strftime("%Y-%m-%d")) |
543
|
|
|
logging.debug("date %s dataset: %s", dateo, sdday_xr) |
544
|
|
|
return sdday_xr |
545
|
|
|
|
546
|
|
|
|
547
|
|
View Code Duplication |
def main(): |
|
|
|
|
548
|
|
|
logging.basicConfig(level=logging.WARNING, |
549
|
|
|
format="[%(levelname)-8s] (%(asctime)s) " |
550
|
|
|
"%(filename)s:%(lineno)d %(message)s", |
551
|
|
|
datefmt="%Y-%m-%d %H:%M:%S %z") |
552
|
|
|
|
553
|
|
|
parser = ap.ArgumentParser() |
554
|
|
|
parser.add_argument("file", default="SCIA_NO.nc", |
555
|
|
|
help="the filename of the output netcdf file") |
556
|
|
|
parser.add_argument("-M", "--month", metavar="YEAR-MM", |
557
|
|
|
help="infer start and end dates for month") |
558
|
|
|
parser.add_argument("-D", "--date_range", metavar="START_DATE:END_DATE", |
559
|
|
|
help="colon-separated start and end dates") |
560
|
|
|
parser.add_argument("-d", "--dates", help="comma-separated list of dates") |
561
|
|
|
parser.add_argument("-B", "--base_date", |
562
|
|
|
metavar="YEAR-MM-DD", default="2000-01-01", |
563
|
|
|
help="Reference date to base the time values (days) on " |
564
|
|
|
"(default: %(default)s).") |
565
|
|
|
parser.add_argument("-f", "--orbit_file", |
566
|
|
|
help="the file containing the input orbits") |
567
|
|
|
parser.add_argument("-r", "--retrieval_version", default="v6.2", |
568
|
|
|
help="SCIAMACHY level 2 data version to process") |
569
|
|
|
parser.add_argument("-R", "--file_version", default="2.3", |
570
|
|
|
help="Postprocessing format version of the output file") |
571
|
|
|
parser.add_argument("-A", "--author", default="unknown", |
572
|
|
|
help="Author of the post-processed data set " |
573
|
|
|
"(default: %(default)s)") |
574
|
|
|
parser.add_argument("-p", "--path", default=None, |
575
|
|
|
help="path containing the L2 data") |
576
|
|
|
parser.add_argument("-s", "--spectra", default=None, metavar="PATH", |
577
|
|
|
help="path containing the L1c spectra") |
578
|
|
|
parser.add_argument("-m", "--mlt", action="store_true", default=False, |
579
|
|
|
help="indicate nominal (False, default) or MLT data (True)") |
580
|
|
|
parser.add_argument("-X", "--xarray", action="store_true", default=False, |
581
|
|
|
help="use xarray to prepare the dataset" |
582
|
|
|
" (experimental, default %(default)s)") |
583
|
|
|
loglevels = parser.add_mutually_exclusive_group() |
584
|
|
|
loglevels.add_argument("-l", "--loglevel", default="WARNING", |
585
|
|
|
choices=['DEBUG', 'INFO', 'WARNING', 'ERROR', 'CRITICAL'], |
586
|
|
|
help="change the loglevel (default: 'WARNING')") |
587
|
|
|
loglevels.add_argument("-q", "--quiet", action="store_true", default=False, |
588
|
|
|
help="less output, same as --loglevel=ERROR (default: False)") |
589
|
|
|
loglevels.add_argument("-v", "--verbose", action="store_true", default=False, |
590
|
|
|
help="verbose output, same as --loglevel=INFO (default: False)") |
591
|
|
|
args = parser.parse_args() |
592
|
|
|
if args.quiet: |
593
|
|
|
logging.getLogger().setLevel(logging.ERROR) |
594
|
|
|
elif args.verbose: |
595
|
|
|
logging.getLogger().setLevel(logging.INFO) |
596
|
|
|
else: |
597
|
|
|
logging.getLogger().setLevel(args.loglevel) |
598
|
|
|
|
599
|
|
|
logging.info("processing L2 version: %s", args.retrieval_version) |
600
|
|
|
logging.info("writing data file version: %s", args.file_version) |
601
|
|
|
|
602
|
|
|
pddrange = [] |
603
|
|
|
if args.month is not None: |
604
|
|
|
d0 = pd.to_datetime(args.month + "-01") |
605
|
|
|
pddrange += pd.date_range(d0, d0 + pd.tseries.offsets.MonthEnd()) |
606
|
|
|
if args.date_range is not None: |
607
|
|
|
pddrange += pd.date_range(*args.date_range.split(':')) |
608
|
|
|
if args.dates is not None: |
609
|
|
|
pddrange += pd.to_datetime(args.dates.split(',')) |
610
|
|
|
logging.debug("pddrange: %s", pddrange) |
611
|
|
|
|
612
|
|
|
olist = [] |
613
|
|
|
for date in pddrange: |
614
|
|
|
try: |
615
|
|
|
olist += get_orbits_from_date(date.strftime("%Y-%m-%d"), |
616
|
|
|
mlt=args.mlt, path=args.path, L2_version=args.retrieval_version) |
617
|
|
|
except: # handle NaT |
618
|
|
|
pass |
619
|
|
|
if args.orbit_file is not None: |
620
|
|
|
olist += np.genfromtxt(args.orbit_file, dtype=np.int32).tolist() |
621
|
|
|
logging.debug("olist: %s", olist) |
622
|
|
|
|
623
|
|
|
if not olist: |
624
|
|
|
logging.warn("No orbits to process.") |
625
|
|
|
return |
626
|
|
|
|
627
|
|
|
sdlist, sdxr_ds = combine_orbit_data(olist, |
628
|
|
|
ref_date=args.base_date, |
629
|
|
|
L2_version=args.retrieval_version, file_version=args.file_version, |
630
|
|
|
dens_path=args.path, spec_base=args.spectra, use_xarray=args.xarray, |
631
|
|
|
save_nc=False) |
632
|
|
|
sdlist.author = args.author |
633
|
|
|
|
634
|
|
|
if args.xarray and sdxr_ds is not None: |
635
|
|
|
sdxr_ds.attrs["author"] = args.author |
636
|
|
|
sd_xr = sddata_xr_set_attrs(sdxr_ds, ref_date=args.base_date) |
637
|
|
|
sd_xr2 = sdlist.to_xarray() |
638
|
|
|
logging.debug(sd_xr) |
639
|
|
|
logging.debug(sd_xr2) |
640
|
|
|
logging.debug("equal datasets: %s", sd_xr.equals(sd_xr2)) |
641
|
|
|
xr.testing.assert_allclose(sd_xr, sd_xr2) |
642
|
|
|
if sd_xr2 is not None: |
643
|
|
|
logging.debug("xarray dataset: %s", sd_xr2) |
644
|
|
|
sd_xr2.to_netcdf(args.file, unlimited_dims=["time"]) |
645
|
|
|
else: |
646
|
|
|
if sdlist.no_dens is not None: |
647
|
|
|
sdlist.write_to_netcdf(args.file) |
648
|
|
|
else: |
649
|
|
|
logging.warn("Processed data is empty.") |
650
|
|
|
|
651
|
|
|
|
652
|
|
|
if __name__ == "__main__": |
653
|
|
|
main() |
654
|
|
|
|