satpy.readers.seviri_l1b_hrit module

SEVIRI HRIT format reader.

Introduction

The seviri_l1b_hrit reader reads and calibrates MSG-SEVIRI L1.5 image data in HRIT format. The format is explained in the MSG Level 1.5 Image Data Format Description. The files are usually named as follows:

H-000-MSG4__-MSG4________-_________-PRO______-201903011200-__
H-000-MSG4__-MSG4________-IR_108___-000001___-201903011200-__
H-000-MSG4__-MSG4________-IR_108___-000002___-201903011200-__
H-000-MSG4__-MSG4________-IR_108___-000003___-201903011200-__
H-000-MSG4__-MSG4________-IR_108___-000004___-201903011200-__
H-000-MSG4__-MSG4________-IR_108___-000005___-201903011200-__
H-000-MSG4__-MSG4________-IR_108___-000006___-201903011200-__
H-000-MSG4__-MSG4________-IR_108___-000007___-201903011200-__
H-000-MSG4__-MSG4________-IR_108___-000008___-201903011200-__
H-000-MSG4__-MSG4________-_________-EPI______-201903011200-__

Each image is decomposed into 24 segments (files) for the high-resolution-visible (HRV) channel and 8 segments for other visible (VIS) and infrared (IR) channels. Additionally, there is one prologue and one epilogue file for the entire scan which contain global metadata valid for all channels.

Reader Arguments

Some arguments can be provided to the reader to change its behaviour. These are provided through the Scene instantiation, eg:

Scene(reader="seviri_l1b_hrit", filenames=fnames, reader_kwargs={'fill_hrv': False})

To see the full list of arguments that can be provided, look into the documentation of HRITMSGFileHandler.

Compression

This reader accepts compressed HRIT files, ending in C_ as other HRIT readers, see satpy.readers.hrit_base.HRITFileHandler.

This reader also accepts bzipped file with the extension .bz2 for the prologue, epilogue, and segment files.

Example

Here is an example how to read the data in satpy:

from satpy import Scene
import glob

filenames = glob.glob('data/H-000-MSG4__-MSG4________-*201903011200*')
scn = Scene(filenames=filenames, reader='seviri_l1b_hrit')
scn.load(['VIS006', 'IR_108'])
print(scn['IR_108'])

Output:

<xarray.DataArray (y: 3712, x: 3712)>
dask.array<shape=(3712, 3712), dtype=float32, chunksize=(464, 3712)>
Coordinates:
    acq_time  (y) datetime64[ns] NaT NaT NaT NaT NaT NaT ... NaT NaT NaT NaT NaT
  * x         (x) float64 5.566e+06 5.563e+06 5.56e+06 ... -5.566e+06 -5.569e+06
  * y         (y) float64 -5.566e+06 -5.563e+06 ... 5.566e+06 5.569e+06
Attributes:
    orbital_parameters:       {'projection_longitude': 0.0, 'projection_latit...
    platform_name:            Meteosat-11
    georef_offset_corrected:  True
    standard_name:            brightness_temperature
    raw_metadata:             {'file_type': 0, 'total_header_length': 6198, '...
    wavelength:               (9.8, 10.8, 11.8)
    units:                    K
    sensor:                   seviri
    platform_name:            Meteosat-11
    start_time:               2019-03-01 12:00:09.716000
    end_time:                 2019-03-01 12:12:42.946000
    area:                     Area ID: some_area_name\\nDescription: On-the-fl...
    name:                     IR_108
    resolution:               3000.403165817
    calibration:              brightness_temperature
    polarization:             None
    level:                    None
    modifiers:                ()
    ancillary_variables:      []

The filenames argument can either be a list of strings, see the example above, or a list of satpy.readers.FSFile objects. FSFiles can be used in conjunction with fsspec, e.g. to handle in-memory data:

import glob

from fsspec.implementations.memory import MemoryFile, MemoryFileSystem
from satpy import Scene
from satpy.readers import FSFile

# In this example, we will make use of `MemoryFile`s in a `MemoryFileSystem`.
memory_fs = MemoryFileSystem()

# Usually, the data already resides in memory.
# For explanatory reasons, we will load the files found with glob in memory,
#  and load the scene with FSFiles.
filenames = glob.glob('data/H-000-MSG4__-MSG4________-*201903011200*')
fs_files = []
for fn in filenames:
    with open(fn, 'rb') as fh:
        fs_files.append(MemoryFile(
            fs=memory_fs,
            path="{}{}".format(memory_fs.root_marker, fn),
            data=fh.read()
        ))
        fs_files[-1].commit()  # commit the file to the filesystem
fs_files = [FSFile(open_file) for open_file in filenames]  # wrap MemoryFiles as FSFiles
# similar to the example above, we pass a list of FSFiles to the `Scene`
scn = Scene(filenames=fs_files, reader='seviri_l1b_hrit')
scn.load(['VIS006', 'IR_108'])
print(scn['IR_108'])

Output:

<xarray.DataArray (y: 3712, x: 3712)>
dask.array<shape=(3712, 3712), dtype=float32, chunksize=(464, 3712)>
Coordinates:
    acq_time  (y) datetime64[ns] NaT NaT NaT NaT NaT NaT ... NaT NaT NaT NaT NaT
  * x         (x) float64 5.566e+06 5.563e+06 5.56e+06 ... -5.566e+06 -5.569e+06
  * y         (y) float64 -5.566e+06 -5.563e+06 ... 5.566e+06 5.569e+06
Attributes:
    orbital_parameters:       {'projection_longitude': 0.0, 'projection_latit...
    platform_name:            Meteosat-11
    georef_offset_corrected:  True
    standard_name:            brightness_temperature
    raw_metadata:             {'file_type': 0, 'total_header_length': 6198, '...
    wavelength:               (9.8, 10.8, 11.8)
    units:                    K
    sensor:                   seviri
    platform_name:            Meteosat-11
    start_time:               2019-03-01 12:00:09.716000
    end_time:                 2019-03-01 12:12:42.946000
    area:                     Area ID: some_area_name\\nDescription: On-the-fl...
    name:                     IR_108
    resolution:               3000.403165817
    calibration:              brightness_temperature
    polarization:             None
    level:                    None
    modifiers:                ()
    ancillary_variables:      []

References

class satpy.readers.seviri_l1b_hrit.HRITMSGEpilogueFileHandler(filename, filename_info, filetype_info, calib_mode='nominal', ext_calib_coefs=None, include_raw_metadata=False, mda_max_array_size=None, fill_hrv=None, mask_bad_quality_scan_lines=None)[source]

Bases: HRITMSGPrologueEpilogueBase

SEVIRI HRIT epilogue reader.

Initialize the reader.

read_epilogue()[source]

Read the epilogue metadata.

reduce(max_size)[source]

Reduce the epilogue metadata.

class satpy.readers.seviri_l1b_hrit.HRITMSGFileHandler(filename, filename_info, filetype_info, prologue, epilogue, calib_mode='nominal', ext_calib_coefs=None, include_raw_metadata=False, mda_max_array_size=100, fill_hrv=True, mask_bad_quality_scan_lines=True)[source]

Bases: HRITFileHandler

SEVIRI HRIT format reader.

Calibration

See satpy.readers.seviri_base.

Padding of the HRV channel

By default, the HRV channel is loaded padded with no-data, that is it is returned as a full-disk dataset. If you want the original, unpadded, data, just provide the fill_hrv as False in the reader_kwargs:

scene = satpy.Scene(filenames,
                    reader='seviri_l1b_hrit',
                    reader_kwargs={'fill_hrv': False})

Metadata

See satpy.readers.seviri_base.

Initialize the reader.

calibrate(data, calibration)[source]

Calibrate the data.

property end_time

Get the end time.

get_area_def(dsid)[source]

Get the area definition of the band.

get_dataset(key, info)[source]

Get the dataset.

property nominal_end_time

Get the end time.

property nominal_start_time

Get the start time.

pad_hrv_data(res)[source]

Add empty pixels around the HRV.

property start_time

Get the start time.

class satpy.readers.seviri_l1b_hrit.HRITMSGPrologueEpilogueBase(filename, filename_info, filetype_info, hdr_info)[source]

Bases: HRITFileHandler

Base reader for prologue and epilogue files.

Initialize the file handler for prologue and epilogue files.

reduce(max_size)[source]

Reduce the metadata (placeholder).

class satpy.readers.seviri_l1b_hrit.HRITMSGPrologueFileHandler(filename, filename_info, filetype_info, calib_mode='nominal', ext_calib_coefs=None, include_raw_metadata=False, mda_max_array_size=None, fill_hrv=None, mask_bad_quality_scan_lines=None)[source]

Bases: HRITMSGPrologueEpilogueBase

SEVIRI HRIT prologue reader.

Initialize the reader.

get_earth_radii()[source]

Get earth radii from prologue.

Returns

Equatorial radius, polar radius [m]

read_prologue()[source]

Read the prologue metadata.

reduce(max_size)[source]

Reduce the prologue metadata.

property satpos

Get actual satellite position in geodetic coordinates (WGS-84).

Evaluate orbit polynomials at the start time of the scan.

Returns: Longitude [deg east], Latitude [deg north] and Altitude [m]

satpy.readers.seviri_l1b_hrit.pad_data(data, final_size, east_bound, west_bound)[source]

Pad the data given east and west bounds and the desired size.