satpy.readers.seviri_base module

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

class satpy.readers.seviri_base.MeirinkCalibrationHandler(calib_mode)[source]

Bases: object

Re-calibration of the SEVIRI visible channels slope (see Meirink 2013).

Initialize the calibration handler.

get_slope(platform, channel, time)[source]

Return the slope using the provided calibration coefficients.

class satpy.readers.seviri_base.MpefProductHeader[source]

Bases: object

MPEF product header class.

get()[source]

Return numpy record_array for MPEF product header.

property images_used

Return structure for images_used.

exception satpy.readers.seviri_base.NoValidOrbitParams[source]

Bases: Exception

Exception when validOrbitParameters are missing.

class satpy.readers.seviri_base.OrbitPolynomial(coefs, start_time, end_time)[source]

Bases: object

Polynomial encoding the satellite position.

Satellite position as a function of time is encoded in the coefficients of an 8th-order Chebyshev polynomial.

Initialize the polynomial.

evaluate(time)[source]

Get satellite position in earth-centered cartesian coordinates.

Parameters:

time – Timestamp where to evaluate the polynomial

Returns:

Earth-centered cartesian coordinates (x, y, z) in meters

class satpy.readers.seviri_base.OrbitPolynomialFinder(orbit_polynomials)[source]

Bases: object

Find orbit polynomial for a given timestamp.

Initialize with the given candidates.

Parameters:

orbit_polynomials

Dictionary of orbit polynomials as found in SEVIRI L1B files:

{'X': x_polynomials,
 'Y': y_polynomials,
 'Z': z_polynomials,
 'StartTime': polynomials_valid_from,
 'EndTime': polynomials_valid_to}

_get_closest_interval(time)[source]

Find interval closest to the given timestamp.

Returns:

Index of closest interval, distance from its center

_get_closest_interval_within(time, threshold)[source]

Find interval closest to the given timestamp within a given distance.

Parameters:
  • time – Timestamp of interest

  • threshold – Maximum distance between timestamp and interval center

Returns:

Index of closest interval

_get_enclosing_interval(time)[source]

Find interval enclosing the given timestamp.

get_orbit_polynomial(time, max_delta=6)[source]

Get orbit polynomial valid for the given time.

Orbit polynomials are only valid for certain time intervals. Find the polynomial, whose corresponding interval encloses the given timestamp. If there are multiple enclosing intervals, use the most recent one. If there is no enclosing interval, find the interval whose centre is closest to the given timestamp (but not more than max_delta hours apart).

Why are there gaps between those intervals? Response from EUM:

A manoeuvre is a discontinuity in the orbit parameters. The flight dynamic algorithms are not made to interpolate over the time-span of the manoeuvre; hence we have elements describing the orbit before a manoeuvre and a new set of elements describing the orbit after the manoeuvre. The flight dynamic products are created so that there is an intentional gap at the time of the manoeuvre. Also the two pre-manoeuvre elements may overlap. But the overlap is not of an issue as both sets of elements describe the same pre-manoeuvre orbit (with negligible variations).

class satpy.readers.seviri_base.SEVIRICalibrationAlgorithm(platform_id, scan_time)[source]

Bases: object

SEVIRI calibration algorithms.

Initialize the calibration algorithm.

_erads2bt(data, channel_name)[source]

Convert effective radiance to brightness temperature.

_srads2bt(data, channel_name)[source]

Convert spectral radiance to brightness temperature.

_tl15(data, wavenumber)[source]

Compute the L15 temperature.

convert_to_radiance(data, gain, offset)[source]

Calibrate to radiance.

ir_calibrate(data, channel_name, cal_type)[source]

Calibrate to brightness temperature.

vis_calibrate(data, solar_irradiance)[source]

Calibrate to reflectance.

This uses the method described in Conversion from radiances to reflectances for SEVIRI warm channels: https://www-cdn.eumetsat.int/files/2020-04/pdf_msg_seviri_rad2refl.pdf

class satpy.readers.seviri_base.SEVIRICalibrationHandler(platform_id, channel_name, coefs, calib_mode, scan_time)[source]

Bases: object

Calibration handler for SEVIRI HRIT-, native- and netCDF-formats.

Handles selection of calibration coefficients and calls the appropriate calibration algorithm.

Initialize the calibration handler.

calibrate(data, calibration)[source]

Calibrate the given data.

get_gain_offset()[source]

Get gain & offset for calibration from counts to radiance.

Choices for internal coefficients are nominal or GSICS. If no GSICS coefficients are available for a certain channel, fall back to nominal coefficients. External coefficients take precedence over internal coefficients.

satpy.readers.seviri_base._create_bad_quality_lines_mask(line_validity, line_geometric_quality, line_radiometric_quality)[source]

Create bad quality scan lines mask.

For details on quality flags see MSG Level 1.5 Image Data Format Description page 109.

Parameters:
  • line_validity (numpy.ndarray) – Quality flags with shape (nlines,).

  • line_geometric_quality (numpy.ndarray) – Quality flags with shape (nlines,).

  • line_radiometric_quality (numpy.ndarray) – Quality flags with shape (nlines,).

Returns:

Indicating if the scan line is bad.

Return type:

numpy.ndarray

satpy.readers.seviri_base.add_scanline_acq_time(dataset, acq_time)[source]

Add scanline acquisition time to the given dataset.

satpy.readers.seviri_base.calculate_area_extent(area_dict)[source]

Calculate the area extent seen by a geostationary satellite.

Parameters:

area_dict – A dictionary containing the required parameters center_point: Center point for the projection north: Northmost row number east: Eastmost column number west: Westmost column number south: Southmost row number column_step: Pixel resolution in meters in east-west direction line_step: Pixel resolution in meters in south-north direction [column_offset: Column offset, defaults to 0 if not given] [line_offset: Line offset, defaults to 0 if not given]

Returns:

An area extent for the scene defined by the lower left and

upper right corners

Return type:

tuple

# For Earth model 2 and full disk VISIR, (center_point - west - 0.5 + we_offset) must be -1856.5 . # See MSG Level 1.5 Image Data Format Description Figure 7 - Alignment and numbering of the non-HRV pixels.

satpy.readers.seviri_base.chebyshev(coefs, time, domain)[source]

Evaluate a Chebyshev Polynomial.

Parameters:
  • coefs (list, np.array) – Coefficients defining the polynomial

  • time (int, float) – Time where to evaluate the polynomial

  • domain (list, tuple) – Domain (or time interval) for which the polynomial is defined: [left, right]

Reference: Appendix A in the MSG Level 1.5 Image Data Format Description.

satpy.readers.seviri_base.chebyshev_3d(coefs, time, domain)[source]

Evaluate Chebyshev Polynomials for three dimensions (x, y, z).

Expects the three coefficient sets to be defined in the same domain.

Parameters:
Returns:

Polynomials evaluated in (x, y, z) dimension.

satpy.readers.seviri_base.create_coef_dict(coefs_nominal, coefs_gsics, radiance_type, ext_coefs)[source]

Create coefficient dictionary expected by calibration class.

satpy.readers.seviri_base.dec10216(inbuf)[source]

Decode 10 bits data into 16 bits words.

/*
 * pack 4 10-bit words in 5 bytes into 4 16-bit words
 *
 * 0       1       2       3       4       5
 * 01234567890123456789012345678901234567890
 * 0         1         2         3         4
 */
ip = &in_buffer[i];
op = &out_buffer[j];
op[0] = ip[0]*4 + ip[1]/64;
op[1] = (ip[1] & 0x3F)*16 + ip[2]/16;
op[2] = (ip[2] & 0x0F)*64 + ip[3]/4;
op[3] = (ip[3] & 0x03)*256 +ip[4];
satpy.readers.seviri_base.get_cds_time(days, msecs)[source]

Compute timestamp given the days since epoch and milliseconds of the day.

1958-01-01 00:00 is interpreted as fill value and will be replaced by NaT (Not a Time).

Parameters:
Returns:

Timestamp(s)

Return type:

numpy.datetime64

satpy.readers.seviri_base.get_meirink_slope(meirink_coefs, acquisition_time)[source]

Compute the slope for the visible channel calibration according to Meirink 2013.

S = A + B * 1.e-3* Day

S is here in µW m-2 sr-1 (cm-1)-1

EUMETSAT calibration is given in mW m-2 sr-1 (cm-1)-1, so an extra factor of 1/1000 must be applied.

satpy.readers.seviri_base.get_padding_area(shape, dtype)[source]

Create a padding area filled with no data.

satpy.readers.seviri_base.get_satpos(orbit_polynomial, time, semi_major_axis, semi_minor_axis)[source]

Get satellite position in geodetic coordinates.

Parameters:
  • orbit_polynomial – OrbitPolynomial instance

  • time – Timestamp where to evaluate the polynomial

  • semi_major_axis – Semi-major axis of the ellipsoid

  • semi_minor_axis – Semi-minor axis of the ellipsoid

Returns:

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

satpy.readers.seviri_base.mask_bad_quality(data, line_validity, line_geometric_quality, line_radiometric_quality)[source]

Mask scan lines with bad quality.

Parameters:
Returns:

data with lines flagged as bad converted to np.nan.

Return type:

xarray.DataArray

satpy.readers.seviri_base.pad_data_horizontally(data, final_size, east_bound, west_bound)[source]

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

satpy.readers.seviri_base.pad_data_vertically(data, final_size, south_bound, north_bound)[source]

Pad the data given south and north bounds and the desired size.

satpy.readers.seviri_base.round_nom_time(dt, time_delta)[source]

Round a datetime object to a multiple of a timedelta.

dt : datetime.datetime object, default now. time_delta : timedelta object, we round to a multiple of this, default 1 minute. adapted for SEVIRI from: https://stackoverflow.com/questions/3463930/how-to-round-the-minute-of-a-datetime-object-python

satpy.readers.seviri_base.should_apply_meirink(calib_mode, channel_name)[source]

Decide whether to use the Meirink calibration coefficients.