Reading

Satpy supports reading and loading data from many input file formats and schemes through the concept of readers. Each reader supports a specific type of input data. The Scene object provides a simple interface around all the complexity of these various formats through its load method. The following sections describe the different way data can be loaded, requested, or added to a Scene object.

Available Readers

For readers currently available in Satpy see Satpy Readers. Additionally to get a list of available readers you can use the available_readers function. By default, it returns the names of available readers. To return additional reader information use available_readers(as_dict=True):

>>> from satpy import available_readers
>>> available_readers()

Filter loaded files

Coming soon…

Load data

Datasets in Satpy are identified by certain pieces of metadata set during data loading. These include name, wavelength, calibration, resolution, polarization, and modifiers. Normally, once a Scene is created requesting datasets by name or wavelength is all that is needed:

>>> from satpy import Scene
>>> scn = Scene(reader="seviri_l1b_hrit", filenames=filenames)
>>> scn.load([0.6, 0.8, 10.8])
>>> scn.load(['IR_120', 'IR_134'])

However, in many cases datasets are available in multiple spatial resolutions, multiple calibrations (brightness_temperature, reflectance, radiance, etc), multiple polarizations, or have corrections or other modifiers already applied to them. By default Satpy will provide the version of the dataset with the highest resolution and the highest level of calibration (brightness temperature or reflectance over radiance). It is also possible to request one of these exact versions of a dataset by using the DataQuery class:

>>> from satpy import DataQuery
>>> my_channel_id = DataQuery(name='IR_016', calibration='radiance')
>>> scn.load([my_channel_id])
>>> print(scn['IR_016'])

Or request multiple datasets at a specific calibration, resolution, or polarization:

>>> scn.load([0.6, 0.8], resolution=1000)

Or multiple calibrations:

>>> scn.load([0.6, 10.8], calibration=['brightness_temperature', 'radiance'])

In the above case Satpy will load whatever dataset is available and matches the specified parameters. So the above load call would load the 0.6 (a visible/reflectance band) radiance data and 10.8 (an IR band) brightness temperature data.

For geostationary satellites that have the individual channel data separated to several files (segments) the missing segments are padded by default to full disk area. This is made to simplify caching of resampling look-up tables (see Resampling for more information). To disable this, the user can pass pad_data keyword argument when loading datasets:

>>> scn.load([0.6, 10.8], pad_data=False)

For geostationary products, where the imagery is stored in the files in an unconventional orientation (e.g. MSG SEVIRI L1.5 data are stored with the southwest corner in the upper right), the keyword argument upper_right_corner can be passed into the load call to automatically flip the datasets to the wished orientation. Accepted argument values are 'NE', 'NW', 'SE', 'SW', and 'native'. By default, no flipping is applied (corresponding to upper_right_corner='native') and the data are delivered in the original format. To get the data in the common upright orientation, load the datasets using e.g.:

>>> scn.load(['VIS008'], upper_right_corner='NE')

Note

If a dataset could not be loaded there is no exception raised. You must check the scn.missing_datasets property for any DataID that could not be loaded.

To find out what datasets are available from a reader from the files that were provided to the Scene use available_dataset_ids():

>>> scn.available_dataset_ids()

Or available_dataset_names() for just the string names of Datasets:

>>> scn.available_dataset_names()

Load remote data

Starting with Satpy version 0.25.1 with supported readers it is possible to load data from remote file systems like s3fs or fsspec. For example:

>>> from satpy import Scene
>>> from satpy.readers import FSFile
>>> import fsspec

>>> filename = 'noaa-goes16/ABI-L1b-RadC/2019/001/17/*_G16_s20190011702186*'

>>> the_files = fsspec.open_files("simplecache::s3://" + filename, s3={'anon': True})

>>> fs_files = [FSFile(open_file) for open_file in the_files]

>>> scn = Scene(filenames=fs_files, reader='abi_l1b')
>>> scn.load(['true_color_raw'])

Check the list of Satpy Readers to see which reader supports remote files. For the usage of fsspec and advanced features like caching files locally see the fsspec Documentation .

Search for local/remote files

Satpy provides a utility find_files_and_readers() for searching for files in a base directory matching various search parameters. This function discovers files based on filename patterns. It returns a dictionary mapping reader name to a list of filenames supported. This dictionary can be passed directly to the Scene initialization.

>>> from satpy import find_files_and_readers, Scene
>>> from datetime import datetime
>>> my_files = find_files_and_readers(base_dir='/data/viirs_sdrs',
...                                   reader='viirs_sdr',
...                                   start_time=datetime(2017, 5, 1, 18, 1, 0),
...                                   end_time=datetime(2017, 5, 1, 18, 30, 0))
>>> scn = Scene(filenames=my_files)

See the find_files_and_readers() documentation for more information on the possible parameters as well as for searching on remote file systems.

Metadata

The datasets held by a scene also provide vital metadata such as dataset name, units, observation time etc. The following attributes are standardized across all readers:

  • name, and other identifying metadata keys: See Satpy internal workings: having a look under the hood.

  • start_time: Left boundary of the time interval covered by the dataset. For more information see the Time Metadata section below.

  • end_time: Right boundary of the time interval covered by the dataset. For more information see the Time Metadata section below.

  • area: AreaDefinition or SwathDefinition if data is geolocated. Areas are used for gridded projected data and Swaths when data must be described by individual longitude/latitude coordinates. See the Coordinates section below.

  • reader: The name of the Satpy reader that produced the dataset.

  • orbital_parameters: Dictionary of orbital parameters describing the satellite’s position. See the Orbital Parameters section below for more information.

  • time_parameters: Dictionary of additional time parameters describing the time ranges related to the requests or schedules for when observations should happen and when they actually do. See Time Metadata below for details.

  • raw_metadata: Raw, unprocessed metadata from the reader.

Note that the above attributes are not necessarily available for each dataset.

Time Metadata

In addition to the generic start_time and end_time pieces of metadata there are other time fields that may be provided if the reader supports them. These items are stored in a time_parameters sub-dictionary and they include values like:

  • observation_start_time: The point in time when a sensor began recording for the current data.

  • observation_end_time: Same as observation_start_time, but when data has stopped being recorded.

  • nominal_start_time: The “human friendly” time describing the start of the data observation interval or repeat cycle. This time is often on a round minute (seconds=0). Along with the nominal end time, these times define the regular interval of the data collection. For example, GOES-16 ABI full disk images are collected every 10 minutes (in the common configuration) so nominal_start_time and nominal_end_time would be 10 minutes apart regardless of when the instrument recorded data inside that interval. This time may also be referred to as the repeat cycle, repeat slot, or time slot.

  • nominal_end_time: Same as nominal_start_time, but the end of the interval.

In general, start_time and end_time will be set to the “nominal” time by the reader. This ensures that other Satpy components get a consistent time for calculations (ex. generation of solar zenith angles) and can be reused between bands.

See the Coordinates section below for more information on time information that may show up as a per-element/row “coordinate” on the DataArray (ex. acquisition time) instead of as metadata.

Orbital Parameters

Orbital parameters describe the position of the satellite. As such they typically come in a few “flavors” for the common types of orbits a satellite may have.

For geostationary satellites it is described using the following scalar attributes:

  • satellite_actual_longitude/latitude/altitude: Current position of the satellite at the time of observation in geodetic coordinates (i.e. altitude is relative and normal to the surface of the ellipsoid). The longitude and latitude are given in degrees, the altitude in meters.

  • satellite_nominal_longitude/latitude/altitude: Center of the station keeping box (a confined area in which the satellite is actively maintained in using maneuvers). Inbetween major maneuvers, when the satellite is permanently moved, the nominal position is constant. The longitude and latitude are given in degrees, the altitude in meters.

  • nadir_longitude/latitude: Intersection of the instrument’s Nadir with the surface of the earth. May differ from the actual satellite position, if the instrument is pointing slightly off the axis (satellite, earth-center). If available, this should be used to compute viewing angles etc. Otherwise, use the actual satellite position. The values are given in degrees.

  • projection_longitude/latitude/altitude: Projection center of the re-projected data. This should be used to compute lat/lon coordinates. Note that the projection center can differ considerably from the actual satellite position. For example MSG-1 was at times positioned at 3.4 degrees west, while the image data was re-projected to 0 degrees. The longitude and latitude are given in degrees, the altitude in meters.

    Note

    For use in pyorbital, the altitude has to be converted to kilometers, see for example pyorbital.orbital.get_observer_look().

For polar orbiting satellites the readers usually provide coordinates and viewing angles of the swath as ancillary datasets. Additional metadata related to the satellite position includes:

  • tle: Two-Line Element (TLE) set used to compute the satellite’s orbit

Coordinates

Each DataArray produced by Satpy has several Xarray coordinate variables added to them.

  • x and y: Projection coordinates for gridded and projected data. By default y and x are the preferred dimensions for all 2D data, but these coordinates are only added for gridded (non-swath) data. For 1D data only the y dimension may be specified.

  • crs: A CRS object defined the Coordinate Reference System for the data. Requires pyproj 2.0 or later to be installed. This is stored as a scalar array by Xarray so it must be accessed by doing crs = my_data_arr.attrs['crs'].item(). For swath data this defaults to a longlat CRS using the WGS84 datum.

  • longitude: Array of longitude coordinates for swath data.

  • latitude: Array of latitude coordinates for swath data.

Readers are free to define any coordinates in addition to the ones above that are automatically added. Other possible coordinates you may see:

  • acq_time: Instrument data acquisition time per scan or row of data.

Adding a Reader to Satpy

This is described in the developer guide, see Adding a Custom Reader to Satpy.

Implemented readers

SEVIRI L1.5 data readers

Common functionality for SEVIRI L1.5 data readers.

Introduction

The Spinning Enhanced Visible and InfraRed Imager (SEVIRI) is the primary instrument on Meteosat Second Generation (MSG) and has the capacity to observe the Earth in 12 spectral channels.

Level 1.5 corresponds to image data that has been corrected for all unwanted radiometric and geometric effects, has been geolocated using a standardised projection, and has been calibrated and radiance-linearised. (From the EUMETSAT documentation)

Satpy provides the following readers for SEVIRI L1.5 data in different formats:

Calibration

This section describes how to control the calibration of SEVIRI L1.5 data.

Calibration to radiance

The SEVIRI L1.5 data readers allow for choosing between two file-internal calibration coefficients to convert counts to radiances:

  • Nominal for all channels (default)

  • GSICS where available (IR currently) and nominal for the remaining channels (VIS & HRV currently)

In order to change the default behaviour, use the reader_kwargs keyword argument upon Scene creation:

import satpy
scene = satpy.Scene(filenames=filenames,
                    reader='seviri_l1b_...',
                    reader_kwargs={'calib_mode': 'GSICS'})
scene.load(['VIS006', 'IR_108'])

In addition, two other calibration methods are available:

  1. It is possible to specify external calibration coefficients for the conversion from counts to radiances. External coefficients take precedence over internal coefficients and over the Meirink coefficients, but you can also mix internal and external coefficients: If external calibration coefficients are specified for only a subset of channels, the remaining channels will be calibrated using the chosen file-internal coefficients (nominal or GSICS). Calibration coefficients must be specified in [mW m-2 sr-1 (cm-1)-1].

  2. The calibration mode meirink-2023 uses coefficients based on an intercalibration with Aqua-MODIS for the visible channels, as found in Inter-calibration of polar imager solar channels using SEVIRI (2013) by J. F. Meirink, R. A. Roebeling, and P. Stammes.

In the following example we use external calibration coefficients for the VIS006 & IR_108 channels, and nominal coefficients for the remaining channels:

coefs = {'VIS006': {'gain': 0.0236, 'offset': -1.20},
         'IR_108': {'gain': 0.2156, 'offset': -10.4}}
scene = satpy.Scene(filenames,
                    reader='seviri_l1b_...',
                    reader_kwargs={'ext_calib_coefs': coefs})
scene.load(['VIS006', 'VIS008', 'IR_108', 'IR_120'])

In the next example we use external calibration coefficients for the VIS006 & IR_108 channels, GSICS coefficients where available (other IR channels) and nominal coefficients for the rest:

coefs = {'VIS006': {'gain': 0.0236, 'offset': -1.20},
         'IR_108': {'gain': 0.2156, 'offset': -10.4}}
scene = satpy.Scene(filenames,
                    reader='seviri_l1b_...',
                    reader_kwargs={'calib_mode': 'GSICS',
                                   'ext_calib_coefs': coefs})
scene.load(['VIS006', 'VIS008', 'IR_108', 'IR_120'])

In the next example we use the mode meirink-2023 calibration coefficients for all visible channels and nominal coefficients for the rest:

scene = satpy.Scene(filenames,
                    reader='seviri_l1b_...',
                    reader_kwargs={'calib_mode': 'meirink-2023'})
scene.load(['VIS006', 'VIS008', 'IR_016'])
Calibration to reflectance

When loading solar channels, the SEVIRI L1.5 data readers apply a correction for the Sun-Earth distance variation throughout the year - as recommended by the EUMETSAT document Conversion from radiances to reflectances for SEVIRI warm channels. In the unlikely situation that this correction is not required, it can be removed on a per-channel basis using satpy.readers.utils.remove_earthsun_distance_correction().

Masking of bad quality scan lines

By default bad quality scan lines are masked and replaced with np.nan for radiance, reflectance and brightness temperature calibrations based on the quality flags provided by the data (for details on quality flags see MSG Level 1.5 Image Data Format Description page 109). To disable masking reader_kwargs={'mask_bad_quality_scan_lines': False} can be passed to the Scene.

Metadata

The SEVIRI L1.5 readers provide the following metadata:

  • The orbital_parameters attribute provides the nominal and actual satellite position, as well as the projection centre. See the Metadata section in the Reading chapter for more information.

  • The acq_time coordinate provides the mean acquisition time for each scanline. Use a MultiIndex to enable selection by acquisition time:

    import pandas as pd
    mi = pd.MultiIndex.from_arrays([scn['IR_108']['y'].data, scn['IR_108']['acq_time'].data],
                                   names=('y_coord', 'time'))
    scn['IR_108']['y'] = mi
    scn['IR_108'].sel(time=np.datetime64('2019-03-01T12:06:13.052000000'))
    
  • Raw metadata from the file header can be included by setting the reader argument include_raw_metadata=True (HRIT and Native format only). Note that this comes with a performance penalty of up to 10% if raw metadata from multiple segments or scans need to be combined. By default, arrays with more than 100 elements are excluded to limit the performance penalty. This threshold can be adjusted using the mda_max_array_size reader keyword argument:

    scene = satpy.Scene(filenames,
                       reader='seviri_l1b_hrit/native',
                       reader_kwargs={'include_raw_metadata': True,
                                      'mda_max_array_size': 1000})
    

References

SEVIRI HRIT format reader

SEVIRI Level 1.5 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:

scn = Scene(filenames=filenames, reader="seviri_l1b_hrit", 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.

Nominal start/end time

Warning

attribute access change

nominal_start_time and nominal_end_time should be accessed using the time_parameters attribute.

nominal_start_time and nominal_end_time are also available directly via start_time and end_time respectively.

Here is an exmaple of the content of the start/end time and time_parameters attibutes

Start time: 2019-08-29 12:00:00
End time:   2019-08-29 12:15:00
time_parameters:
                {'nominal_start_time': datetime.datetime(2019, 8, 29, 12, 0),
                 'nominal_end_time': datetime.datetime(2019, 8, 29, 12, 15),
                 'observation_start_time': datetime.datetime(2019, 8, 29, 12, 0, 9, 338000),
                 'observation_end_time': datetime.datetime(2019, 8, 29, 12, 15, 9, 203000)
                 }
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

SEVIRI Native format reader

SEVIRI Level 1.5 native format reader.

Introduction

The seviri_l1b_native reader reads and calibrates MSG-SEVIRI L1.5 image data in binary format. The format is explained in the MSG Level 1.5 Native Format File Definition. The files are usually named as follows:

MSG4-SEVI-MSG15-0100-NA-20210302124244.185000000Z-NA.nat
Reader Arguments

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

scn = Scene(filenames=filenames, reader="seviri_l1b_native", reader_kwargs={'fill_disk': True})

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

Example:

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

NOTE: When loading the data, the orientation of the image can be set with upper_right_corner-keyword. Possible options are NW, NE, SW, SE, or native.

from satpy import Scene

filenames = ['MSG4-SEVI-MSG15-0100-NA-20210302124244.185000000Z-NA.nat']
scn = Scene(filenames=filenames, reader='seviri_l1b_native')
scn.load(['VIS006', 'IR_108'], upper_right_corner='NE')
print(scn['IR_108'])

Output:

<xarray.DataArray 'reshape-969ef97d34b7b0c70ca19f53c6abcb68' (y: 3712, x: 3712)>
dask.array<truediv, shape=(3712, 3712), dtype=float32, chunksize=(928, 3712), chunktype=numpy.ndarray>
Coordinates:
    acq_time  (y) datetime64[ns] NaT NaT NaT NaT NaT NaT ... NaT NaT NaT NaT NaT
    crs       object PROJCRS["unknown",BASEGEOGCRS["unknown",DATUM["unknown",...
  * y         (y) float64 -5.566e+06 -5.563e+06 ... 5.566e+06 5.569e+06
  * x         (x) float64 5.566e+06 5.563e+06 5.56e+06 ... -5.566e+06 -5.569e+06
Attributes:
    orbital_parameters:       {'projection_longitude': 0.0, 'projection_latit...
    time_parameters:          {'nominal_start_time': datetime.datetime(2021, ...
    units:                    K
    wavelength:               10.8 µm (9.8-11.8 µm)
    standard_name:            toa_brightness_temperature
    platform_name:            Meteosat-11
    sensor:                   seviri
    georef_offset_corrected:  True
    start_time:               2021-03-02 12:30:11.584603
    end_time:                 2021-03-02 12:45:09.949762
    reader:                   seviri_l1b_native
    area:                     Area ID: msg_seviri_fes_3km\\nDescription: MSG S...
    name:                     IR_108
    resolution:               3000.403165817
    calibration:              brightness_temperature
    modifiers:                ()
    _satpy_id:                DataID(name='IR_108', wavelength=WavelengthRang...
    ancillary_variables:      []

References

SEVIRI netCDF format reader

SEVIRI netcdf format reader.

Other xRIT-based readers

HRIT/LRIT format reader.

This module is the base module for all HRIT-based formats. Here, you will find the common building blocks for hrit reading.

One of the features here is the on-the-fly decompression of hrit files. It needs a path to the xRITDecompress binary to be provided through the environment variable called XRIT_DECOMPRESS_PATH. When compressed hrit files are then encountered (files finishing with .C_), they are decompressed to the system’s temporary directory for reading.

JMA HRIT format reader

HRIT format reader for JMA data.

Introduction

The JMA HRIT format is described in the JMA HRIT - Mission Specific Implementation. There are three readers for this format in Satpy:

  • jami_hrit: For data from the JAMI instrument on MTSAT-1R

  • mtsat2-imager_hrit: For data from the Imager instrument on MTSAT-2

  • ahi_hrit: For data from the AHI instrument on Himawari-8/9

Although the data format is identical, the instruments have different characteristics, which is why there is a dedicated reader for each of them. Sample data is available here:

Example:

Here is an example how to read Himwari-8 HRIT data with Satpy:

from satpy import Scene
import glob

filenames = glob.glob('data/IMG_DK01B14_2018011109*')
scn = Scene(filenames=filenames, reader='ahi_hrit')
scn.load(['B14'])
print(scn['B14'])

Output:

<xarray.DataArray (y: 5500, x: 5500)>
dask.array<concatenate, shape=(5500, 5500), dtype=float64, chunksize=(550, 4096), ...
Coordinates:
    acq_time  (y) datetime64[ns] 2018-01-11T09:00:20.995200 ... 2018-01-11T09:09:40.348800
    crs       object +proj=geos +lon_0=140.7 +h=35785831 +x_0=0 +y_0=0 +a=6378169 ...
  * y         (y) float64 5.5e+06 5.498e+06 5.496e+06 ... -5.496e+06 -5.498e+06
  * x         (x) float64 -5.498e+06 -5.496e+06 -5.494e+06 ... 5.498e+06 5.5e+06
Attributes:
    orbital_parameters:   {'projection_longitude': 140.7, 'projection_latitud...
    standard_name:        toa_brightness_temperature
    level:                None
    wavelength:           (11.0, 11.2, 11.4)
    units:                K
    calibration:          brightness_temperature
    file_type:            ['hrit_b14_seg', 'hrit_b14_fd']
    modifiers:            ()
    polarization:         None
    sensor:               ahi
    name:                 B14
    platform_name:        Himawari-8
    resolution:           4000
    start_time:           2018-01-11 09:00:20.995200
    end_time:             2018-01-11 09:09:40.348800
    area:                 Area ID: FLDK, Description: Full Disk, Projection I...
    ancillary_variables:  []

JMA HRIT data contain the scanline acquisition time for only a subset of scanlines. Timestamps of the remaining scanlines are computed using linear interpolation. This is what you’ll find in the acq_time coordinate of the dataset.

Compression

Gzip-compressed MTSAT files can be decompressed on the fly using FSFile:

import fsspec
from satpy import Scene
from satpy.readers import FSFile

filename = "/data/HRIT_MTSAT1_20090101_0630_DK01IR1.gz"
open_file = fsspec.open(filename, compression="gzip")
fs_file = FSFile(open_file)
scn = Scene([fs_file], reader="jami_hrit")
scn.load(["IR1"])

GOES HRIT format reader

GOES HRIT format reader.

References

LRIT/HRIT Mission Specific Implementation, February 2012 GVARRDL98.pdf 05057_SPE_MSG_LRIT_HRI

Electro-L HRIT format reader

HRIT format reader.

References

ELECTRO-L GROUND SEGMENT MSU-GS INSTRUMENT,

LRIT/HRIT Mission Specific Implementation, February 2012

hdf-eos based readers

Modis level 1b hdf-eos format reader.

Introduction

The modis_l1b reader reads and calibrates Modis L1 image data in hdf-eos format. Files often have a pattern similar to the following one:

M[O/Y]D02[1/H/Q]KM.A[date].[time].[collection].[processing_time].hdf

Other patterns where “collection” and/or “proccessing_time” are missing might also work (see the readers yaml file for details). Geolocation files (MOD03) are also supported. The IMAPP direct broadcast naming format is also supported with names like: a1.12226.1846.1000m.hdf.

Saturation Handling

Band 2 of the MODIS sensor is available in 250m, 500m, and 1km resolutions. The band data may include a special fill value to indicate when the detector was saturated in the 250m version of the data. When the data is aggregated to coarser resolutions this saturation fill value is converted to a “can’t aggregate” fill value. By default, Satpy will replace these fill values with NaN to indicate they are invalid. This is typically undesired when generating images for the data as they appear as “holes” in bright clouds. To control this the keyword argument mask_saturated can be passed and set to False to set these two fill values to the maximum valid value.

scene = satpy.Scene(filenames=filenames,
                    reader='modis_l1b',
                    reader_kwargs={'mask_saturated': False})
scene.load(['2'])

Note that the saturation fill value can appear in other bands (ex. bands 7-19) in addition to band 2. Also, the “can’t aggregate” fill value is a generic “catch all” for any problems encountered when aggregating high resolution bands to lower resolutions. Filling this with the max valid value could replace non-saturated invalid pixels with valid values.

Geolocation files

For the 1km data (mod021km) geolocation files (mod03) are optional. If not given to the reader 1km geolocations will be interpolated from the 5km geolocation contained within the file.

For the 500m and 250m data geolocation files are needed.

References

Modis level 2 hdf-eos format reader.

Introduction

The modis_l2 reader reads and calibrates Modis L2 image data in hdf-eos format. Since there are a multitude of different level 2 datasets not all of theses are implemented (yet).

Currently the reader supports:
  • m[o/y]d35_l2: cloud_mask dataset

  • some datasets in m[o/y]d06 files

To get a list of the available datasets for a given file refer to the “Load data” section in Reading.

Geolocation files

Similar to the modis_l1b reader the geolocation files (mod03) for the 1km data are optional and if not given 1km geolocations will be interpolated from the 5km geolocation contained within the file.

For the 500m and 250m data geolocation files are needed.

References

satpy cf nc readers

Reader for files produced with the cf netcdf writer in satpy.

Introduction

The satpy_cf_nc reader reads data written by the satpy cf_writer. Filenames for cf_writer are optional. There are several readers using the same satpy_cf_nc.py reader.

  • Generic reader satpy_cf_nc

  • EUMETSAT GAC FDR reader avhrr_l1c_eum_gac_fdr_nc

Generic reader

The generic satpy_cf_nc reader reads files of type:

'{platform_name}-{sensor}-{start_time:%Y%m%d%H%M%S}-{end_time:%Y%m%d%H%M%S}.nc'

Example:

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

from satpy import Scene

filenames = ['data/npp-viirs-mband-20201007075915-20201007080744.nc']
scn = Scene(reader='satpy_cf_nc', filenames=filenames)
scn.load(['M05'])
scn['M05']

Output:

<xarray.DataArray 'M05' (y: 4592, x: 3200)>
dask.array<open_dataset-d91cfbf1bf4f14710d27446d91cdc6e4M05, shape=(4592, 3200),
    dtype=float32, chunksize=(4096, 3200), chunktype=numpy.ndarray>
Coordinates:
    longitude  (y, x) float32 dask.array<chunksize=(4096, 3200), meta=np.ndarray>
    latitude   (y, x) float32 dask.array<chunksize=(4096, 3200), meta=np.ndarray>
Dimensions without coordinates: y, x
Attributes:
    start_time:                   2020-10-07 07:59:15
    start_orbit:                  46350
    end_time:                     2020-10-07 08:07:44
    end_orbit:                    46350
    calibration:                  reflectance
    long_name:                    M05
    modifiers:                    ('sunz_corrected',)
    platform_name:                Suomi-NPP
    resolution:                   742
    sensor:                       viirs
    standard_name:                toa_bidirectional_reflectance
    units:                        %
    wavelength:                   0.672 µm (0.662-0.682 µm)
    date_created:                 2020-10-07T08:20:02Z
    instrument:                   VIIRS

Notes

Available datasets and attributes will depend on the data saved with the cf_writer.

EUMETSAT AVHRR GAC FDR L1C reader

The avhrr_l1c_eum_gac_fdr_nc reader reads files of type:

''AVHRR-GAC_FDR_1C_{platform}_{start_time:%Y%m%dT%H%M%SZ}_{end_time:%Y%m%dT%H%M%SZ}_{processing_mode}_{disposition_mode}_{creation_time}_{version_int:04d}.nc'

Example:

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

from satpy import Scene

filenames = ['data/AVHRR-GAC_FDR_1C_N06_19810330T042358Z_19810330T060903Z_R_O_20200101T000000Z_0100.nc']
scn = Scene(reader='avhrr_l1c_eum_gac_fdr_nc', filenames=filenames)
scn.load(['brightness_temperature_channel_4'])
scn['brightness_temperature_channel_4']

Output:

<xarray.DataArray 'brightness_temperature_channel_4' (y: 11, x: 409)>
dask.array<open_dataset-55ffbf3623b32077c67897f4283640a5brightness_temperature_channel_4, shape=(11, 409),
    dtype=float32, chunksize=(11, 409), chunktype=numpy.ndarray>
Coordinates:
  * x          (x) int16 0 1 2 3 4 5 6 7 8 ... 401 402 403 404 405 406 407 408
  * y          (y) int64 0 1 2 3 4 5 6 7 8 9 10
    acq_time   (y) datetime64[ns] dask.array<chunksize=(11,), meta=np.ndarray>
    longitude  (y, x) float64 dask.array<chunksize=(11, 409), meta=np.ndarray>
    latitude   (y, x) float64 dask.array<chunksize=(11, 409), meta=np.ndarray>
Attributes:
    start_time:                            1981-03-30 04:23:58
    end_time:                              1981-03-30 06:09:03
    calibration:                           brightness_temperature
    modifiers:                             ()
    resolution:                            1050
    standard_name:                         toa_brightness_temperature
    units:                                 K
    wavelength:                            10.8 µm (10.3-11.3 µm)
    Conventions:                           CF-1.8 ACDD-1.3
    comment:                               Developed in cooperation with EUME...
    creator_email:                         ops@eumetsat.int
    creator_name:                          EUMETSAT
    creator_url:                           https://www.eumetsat.int/
    date_created:                          2020-09-14T10:50:51.073707
    disposition_mode:                      O
    gac_filename:                          NSS.GHRR.NA.D81089.S0423.E0609.B09...
    geospatial_lat_max:                    89.95386902434623
    geospatial_lat_min:                    -89.97581969005503
    geospatial_lat_resolution:             1050 meters
    geospatial_lat_units:                  degrees_north
    geospatial_lon_max:                    179.99952992568998
    geospatial_lon_min:                    -180.0
    geospatial_lon_resolution:             1050 meters
    geospatial_lon_units:                  degrees_east
    ground_station:                        GC
    id:                                    DOI:10.5676/EUM/AVHRR_GAC_L1C_FDR/...
    institution:                           EUMETSAT
    instrument:                            Earth Remote Sensing Instruments >...
    keywords:                              ATMOSPHERE > ATMOSPHERIC RADIATION...
    keywords_vocabulary:                   GCMD Science Keywords, Version 9.1
    licence:                               EUMETSAT data policy https://www.e...
    naming_authority:                      int.eumetsat
    orbit_number_end:                      9123
    orbit_number_start:                    9122
    orbital_parameters_tle:                ['1 11416U 79057A   81090.16350942...
    platform:                              Earth Observation Satellites > NOA...
    processing_level:                      1C
    processing_mode:                       R
    product_version:                       1.0.0
    references:                            Devasthale, A., M. Raspaud, C. Sch...
    source:                                AVHRR GAC Level 1 Data
    standard_name_vocabulary:              CF Standard Name Table v73
    summary:                               Fundamental Data Record (FDR) of m...
    sun_earth_distance_correction_factor:  0.9975244779999585
    time_coverage_end:                     19820803T003900Z
    time_coverage_start:                   19800101T000000Z
    title:                                 AVHRR GAC L1C FDR
    version_calib_coeffs:                  PATMOS-x, v2017r1
    version_pygac:                         1.4.0
    version_pygac_fdr:                     0.1.dev107+gceb7b26.d20200910
    version_satpy:                         0.21.1.dev894+g5cf76e6
    history:                               Created by pytroll/satpy on 2020-0...
    name:                                  brightness_temperature_channel_4
    _satpy_id:                             DataID(name='brightness_temperatur...
    ancillary_variables:                   []

hdf5 based readers

Advanced Geostationary Radiation Imager reader for the Level_1 HDF format.

The files read by this reader are described in the official Real Time Data Service:

Geostationary High-speed Imager reader for the Level_1 HDF format.

This instrument is aboard the Fengyun-4B satellite. No document is available to describe this format is available, but it’s broadly similar to the co-flying AGRI instrument.

Arctica-M N1 HDF5 format reader

Reader for the Arctica-M1 MSU-GS/A data.

The files for this reader are HDF5 and contain channel data at 1km resolution for the VIS channels and 4km resolution for the IR channels. Geolocation data is available at both resolutions, as is sun and satellite geometry.

This reader was tested on sample data provided by EUMETSAT.