Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add meta_data dictionary to reader class. #41

Merged
Merged
Show file tree
Hide file tree
Changes from 4 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
8 changes: 8 additions & 0 deletions pygac/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -22,6 +22,7 @@

import logging
import os
import numpy as np

from pygac.version import __version__ # noqa

Expand Down Expand Up @@ -52,3 +53,10 @@ def centered_modulus(array, divisor):
arr = array % divisor
arr[arr > divisor / 2] -= divisor
return arr


def calculate_sun_earth_distance_correction(jday):
"""Calculate the sun earth distance correction."""
sfinkens marked this conversation as resolved.
Show resolved Hide resolved
# Earth-Sun distance correction factor
corr = 1.0 - 0.0334 * np.cos(2.0 * np.pi * (jday - 2) / 365.25)
return corr
10 changes: 5 additions & 5 deletions pygac/gac_io.py
Original file line number Diff line number Diff line change
Expand Up @@ -83,7 +83,11 @@ def save_gac(satellite_name,
bt3, bt4, bt5,
sun_zen, sat_zen, sun_azi, sat_azi, rel_azi,
qual_flags, start_line, end_line,
gac_file, midnight_scanline, miss_lines):
gac_file, meta_data):

midnight_scanline = meta_data['midnight_scanline']
miss_lines = meta_data['missing_scanlines']
corr = meta_data['sun_earth_distance_correction_factor']

last_scan_line_number = qual_flags[-1, 0]

Expand Down Expand Up @@ -199,10 +203,6 @@ def save_gac(satellite_name,
starttime = start.strftime("%H%M%S%f")[:-5]
enddate = end.strftime("%Y%m%d")
endtime = end.strftime("%H%M%S%f")[:-5]
jday = int(start.strftime("%j"))

# Earth-Sun distance correction factor
corr = 1.0 - 0.0334 * np.cos(2.0 * np.pi * (jday - 2) / 365.25)

# Apply scaling & offset
bt3 -= 273.15
Expand Down
3 changes: 1 addition & 2 deletions pygac/gac_klm.py
Original file line number Diff line number Diff line change
Expand Up @@ -624,8 +624,7 @@ def main(filename, start_line, end_line):
sun_zen, sat_zen, sun_azi, sat_azi, rel_azi,
qual_flags, start_line, end_line,
reader.filename,
reader.get_midnight_scanline(),
reader.get_miss_lines())
reader.meta_data)
LOG.info("pygac took: %s", str(datetime.datetime.now() - tic))


Expand Down
3 changes: 1 addition & 2 deletions pygac/gac_pod.py
Original file line number Diff line number Diff line change
Expand Up @@ -430,8 +430,7 @@ def main(filename, start_line, end_line):
sun_zen, sat_zen, sun_azi, sat_azi, rel_azi,
qual_flags, start_line, end_line,
reader.filename,
reader.get_midnight_scanline(),
reader.get_miss_lines())
reader.meta_data)
LOG.info("pygac took: %s", str(datetime.datetime.now() - tic))


Expand Down
26 changes: 24 additions & 2 deletions pygac/gac_reader.py
Original file line number Diff line number Diff line change
Expand Up @@ -31,6 +31,7 @@
import types

from pygac import (CONFIG_FILE, centered_modulus,
calculate_sun_earth_distance_correction,
get_absolute_azimuth_angle_diff)
try:
import ConfigParser
Expand Down Expand Up @@ -64,6 +65,7 @@ def __init__(self, interpolate_coords=True, adjust_clock_drift=True,
tle_thresh: Maximum number of days between observation and nearest
TLE
"""
self.meta_data = {}
self.interpolate_coords = interpolate_coords
self.adjust_clock_drift = adjust_clock_drift
self.tle_dir = tle_dir
Expand Down Expand Up @@ -255,15 +257,35 @@ def compute_lonlat(self, width, utcs=None, clock_drift_adjust=True):

return lons.reshape(-1, width), lats.reshape(-1, width)

def get_sun_earth_distance_correction(self):
"""Get the julian day and the sun-earth distance correction."""
self.get_times()
year = self.times[0].year
delta = self.times[0].date() - datetime.date(year, 1, 1)
jday = delta.days + 1
sfinkens marked this conversation as resolved.
Show resolved Hide resolved
return calculate_sun_earth_distance_correction(jday)

def update_meta_data(self):
"""Add some metd data to the meta_data dicitonary."""
if 'sun_earth_distance_correction_factor' not in self.meta_data.keys():
sfinkens marked this conversation as resolved.
Show resolved Hide resolved
self.meta_data['sun_earth_distance_correction_factor'] = (
self.get_sun_earth_distance_correction())
if 'midnight_scanline' not in self.meta_data.keys():
self.meta_data['midnight_scanline'] = self.get_midnight_scanline()
if 'missing_scanlines' not in self.meta_data.keys():
self.meta_data['missing_scanlines'] = (
sfinkens marked this conversation as resolved.
Show resolved Hide resolved
self.get_miss_lines()).astype(int)

def get_calibrated_channels(self):
"""Calibrate the solar channels."""
channels = self.get_counts()
self.get_times()
self.update_meta_data()
year = self.times[0].year
delta = self.times[0].date() - datetime.date(year, 1, 1)
jday = delta.days + 1

# Earth-Sun distance correction factor
corr = 1.0 - 0.0334 * np.cos(2.0 * np.pi * (jday - 2) / 365.25)
corr = self.meta_data['sun_earth_distance_correction_factor']

# how many reflective channels are there ?
tot_ref = channels.shape[2] - 3
Expand Down
3 changes: 1 addition & 2 deletions pygac/tests/test_io.py
Original file line number Diff line number Diff line change
Expand Up @@ -212,8 +212,7 @@ def test_save_gac(self, check_user_scanlines, slice_channel, avhrr_gac_io,
rel_azi=mm,
qual_flags=mm,
gac_file=mm,
midnight_scanline=mm,
miss_lines=mm
meta_data=mm
)
slice_channel.return_value = mm, 'miss', 'midnight'
strip_invalid_lat.return_value = 0, 0
Expand Down
9 changes: 9 additions & 0 deletions pygac/tests/test_reader.py
Original file line number Diff line number Diff line change
Expand Up @@ -341,6 +341,15 @@ def test_correct_times_thresh(self, get_header_timestamp):
self.reader.correct_times_thresh()
numpy.testing.assert_array_equal(self.reader.utcs, utcs_expected)

def test_calculate_sun_earth_distance_correction(self):
"""Test the calculate sun earth distance correction method."""
self.reader.utcs = np.array([315748035469, 315748359969,
315751135469, 315754371969,
315754371969]).astype('datetime64[ms]')
self.reader.times = self.reader.to_datetime(self.reader.utcs)
corr = self.reader.get_sun_earth_distance_correction()
numpy.testing.assert_almost_equal(corr, 0.96660494, decimal=7)


def suite():
"""Test suite for test_reader."""
Expand Down