Python interface to map GRIB files to the NetCDF Common Data Model following the CF Convention using ecCodes.
Project description
cfgrib: A Python interface to map GRIB files to the NetCDF Common Data Model following the CF Convention using ecCodes
Python interface to map GRIB files to the Unidata’s Common Data Model v4 following the CF Conventions. The high level API is designed to support a GRIB engine for xarray and it is inspired by netCDF4-python and h5netcdf. Low level access and decoding is performed via the ECMWF ecCodes library and the eccodes python package.
Features with development status Beta:
enables the engine='cfgrib' option to read GRIB files with xarray,
reads most GRIB 1 and 2 files including heterogeneous ones with cfgrib.open_datasets,
supports all modern versions of Python 3.9, 3.8, 3.7 and PyPy3,
the 0.9.6.x series with support for Python 2 will stay active and receive critical bugfixes,
works wherever eccodes-python does: Linux, MacOS and Windows
conda-forge package on all supported platforms,
reads the data lazily and efficiently in terms of both memory usage and disk access,
allows larger-than-memory and distributed processing via xarray and dask,
supports translating coordinates to different data models and naming conventions,
supports writing the index of a GRIB file to disk, to save a full-file scan on open,
accepts objects implementing a generic Fieldset interface as described in ADVANCED_USAGE.rst.
Work in progress:
Beta install a cfgrib utility that can convert a GRIB file to_netcdf with a optional conversion to a specific coordinates data model, see #40.
Alpha/Broken support writing carefully-crafted xarray.Dataset’s to a GRIB1 or GRIB2 file, see the Advanced write usage section below, #18 and #156.
Limitations:
relies on ecCodes for the CF attributes of the data variables,
relies on ecCodes for anything related to coordinate systems / gridType, see #28.
Installation
The easiest way to install cfgrib and all its binary dependencies is via Conda:
$ conda install -c conda-forge cfgrib
alternatively, if you install the binary dependencies yourself, you can install the Python package from PyPI with:
$ pip install cfgrib
Binary dependencies
cfgrib depends on the eccodes python package to access the ECMWF ecCodes binary library, when not using conda please follow the System dependencies section there.
You may run a simple selfcheck command to ensure that your system is set up correctly:
$ python -m cfgrib selfcheck Found: ecCodes v2.20.0. Your system is ready.
Usage
First, you need a well-formed GRIB file, if you don’t have one at hand you can download our ERA5 on pressure levels sample:
$ wget http://download.ecmwf.int/test-data/cfgrib/era5-levels-members.grib
Read-only xarray GRIB engine
Most of cfgrib users want to open a GRIB file as a xarray.Dataset and need to have xarray installed:
$ pip install xarray
In a Python interpreter try:
>>> import xarray as xr >>> ds = xr.open_dataset('era5-levels-members.grib', engine='cfgrib') >>> ds <xarray.Dataset> Dimensions: (number: 10, time: 4, isobaricInhPa: 2, latitude: 61, longitude: 120) Coordinates: * number (number) int64 0 1 2 3 4 5 6 7 8 9 * time (time) datetime64[ns] 2017-01-01 ... 2017-01-02T12:00:00 step timedelta64[ns] ... * isobaricInhPa (isobaricInhPa) float64 850.0 500.0 * latitude (latitude) float64 90.0 87.0 84.0 81.0 ... -84.0 -87.0 -90.0 * longitude (longitude) float64 0.0 3.0 6.0 9.0 ... 351.0 354.0 357.0 valid_time (time) datetime64[ns] ... Data variables: z (number, time, isobaricInhPa, latitude, longitude) float32 ... t (number, time, isobaricInhPa, latitude, longitude) float32 ... Attributes: GRIB_edition: 1 GRIB_centre: ecmf GRIB_centreDescription: European Centre for Medium-Range Weather Forecasts GRIB_subCentre: 0 Conventions: CF-1.7 institution: European Centre for Medium-Range Weather Forecasts history: ...
The cfgrib engine supports all read-only features of xarray like:
merge the content of several GRIB files into a single dataset using xarray.open_mfdataset,
work with larger-than-memory datasets with dask,
allow distributed processing with dask.distributed.
Read arbitrary GRIB keys
By default cfgrib reads a limited set of ecCodes recognised keys from the GRIB files and exposes them as Dataset or DataArray attributes with the GRIB_ prefix. It is possible to have cfgrib read additional keys to the attributes by adding the read_keys dictionary key to the backend_kwargs with values the list of desired GRIB keys:
>>> ds = xr.open_dataset('era5-levels-members.grib', engine='cfgrib', ... backend_kwargs={'read_keys': ['experimentVersionNumber']}) >>> ds.t.attrs['GRIB_experimentVersionNumber'] '0001'
Translate to a custom data model
Contrary to netCDF the GRIB data format is not self-describing and several details of the mapping to the Unidata Common Data Model are arbitrarily set by the software components decoding the format. Details like names and units of the coordinates are particularly important because xarray broadcast and selection rules depend on them. cf2cfm is a small coordinate translation module distributed with cfgrib that make it easy to translate CF compliant coordinates, like the one provided by cfgrib, to a user-defined custom data model with set out_name, units and stored_direction.
For example to translate a cfgrib styled xr.Dataset to the classic ECMWF coordinate naming conventions you can:
>>> import cf2cdm >>> ds = xr.open_dataset('era5-levels-members.grib', engine='cfgrib') >>> cf2cdm.translate_coords(ds, cf2cdm.ECMWF) <xarray.Dataset> Dimensions: (number: 10, time: 4, level: 2, latitude: 61, longitude: 120) Coordinates: * number (number) int64 0 1 2 3 4 5 6 7 8 9 * time (time) datetime64[ns] 2017-01-01 ... 2017-01-02T12:00:00 step timedelta64[ns] ... * level (level) float64 850.0 500.0 * latitude (latitude) float64 90.0 87.0 84.0 81.0 ... -84.0 -87.0 -90.0 * longitude (longitude) float64 0.0 3.0 6.0 9.0 ... 348.0 351.0 354.0 357.0 valid_time (time) datetime64[ns] ... Data variables: z (number, time, level, latitude, longitude) float32 ... t (number, time, level, latitude, longitude) float32 ... Attributes: GRIB_edition: 1 GRIB_centre: ecmf GRIB_centreDescription: European Centre for Medium-Range Weather Forecasts GRIB_subCentre: 0 Conventions: CF-1.7 institution: European Centre for Medium-Range Weather Forecasts history: ...
To translate to the Common Data Model of the Climate Data Store use:
>>> import cf2cdm >>> cf2cdm.translate_coords(ds, cf2cdm.CDS) <xarray.Dataset> Dimensions: (realization: 10, forecast_reference_time: 4, plev: 2, lat: 61, lon: 120) Coordinates: * realization (realization) int64 0 1 2 3 4 5 6 7 8 9 * forecast_reference_time (forecast_reference_time) datetime64[ns] 2017-01... leadtime timedelta64[ns] ... * plev (plev) float64 8.5e+04 5e+04 * lat (lat) float64 -90.0 -87.0 -84.0 ... 84.0 87.0 90.0 * lon (lon) float64 0.0 3.0 6.0 9.0 ... 351.0 354.0 357.0 time (forecast_reference_time) datetime64[ns] ... Data variables: z (realization, forecast_reference_time, plev, lat, lon) float32 ... t (realization, forecast_reference_time, plev, lat, lon) float32 ... Attributes: GRIB_edition: 1 GRIB_centre: ecmf GRIB_centreDescription: European Centre for Medium-Range Weather Forecasts GRIB_subCentre: 0 Conventions: CF-1.7 institution: European Centre for Medium-Range Weather Forecasts history: ...
Filter heterogeneous GRIB files
xr.open_dataset can open a GRIB file only if all the messages with the same shortName can be represented as a single hypercube. For example, a variable t cannot have both isobaricInhPa and hybrid typeOfLevel’s, as this would result in multiple hypercubes for the same variable. Opening a non-conformant GRIB file will fail with a ValueError: multiple values for unique key... error message, see #2.
Furthermore if different variables depend on the same coordinate, for example step, the values of the coordinate must match exactly. For example, if variables t and z share the same step coordinate, they must both have exactly the same set of steps. Opening a non-conformant GRIB file will fail with a ValueError: key present and new value is different... error message, see #13.
In most cases you can handle complex GRIB files containing heterogeneous messages by passing the filter_by_keys key in backend_kwargs to select which GRIB messages belong to a well formed set of hypercubes.
For example to open US National Weather Service complex GRIB2 files you can use:
>>> xr.open_dataset('nam.t00z.awp21100.tm00.grib2', engine='cfgrib', ... backend_kwargs={'filter_by_keys': {'typeOfLevel': 'surface'}}) <xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] ... step timedelta64[ns] ... surface float64 ... latitude (y, x) float64 ... longitude (y, x) float64 ... valid_time datetime64[ns] ... Dimensions without coordinates: y, x Data variables: gust (y, x) float32 ... sp (y, x) float32 ... orog (y, x) float32 ... tp (y, x) float32 ... acpcp (y, x) float32 ... csnow (y, x) float32 ... cicep (y, x) float32 ... cfrzr (y, x) float32 ... crain (y, x) float32 ... cape (y, x) float32 ... cin (y, x) float32 ... hpbl (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP... history: ... >>> xr.open_dataset('nam.t00z.awp21100.tm00.grib2', engine='cfgrib', ... backend_kwargs={'filter_by_keys': {'typeOfLevel': 'heightAboveGround', 'level': 2}}) <xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] ... step timedelta64[ns] ... heightAboveGround float64 ... latitude (y, x) float64 ... longitude (y, x) float64 ... valid_time datetime64[ns] ... Dimensions without coordinates: y, x Data variables: t2m (y, x) float32 ... r2 (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP... history: ...
Automatic filtering
cfgrib also provides a function that automates the selection of appropriate filter_by_keys and returns a list of all valid xarray.Dataset’s in the GRIB file.
>>> import cfgrib >>> cfgrib.open_datasets('nam.t00z.awp21100.tm00.grib2') [<xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 atmosphereSingleLayer float64 0.0 latitude (y, x) float64 ... longitude (y, x) float64 ... valid_time datetime64[ns] ... Dimensions without coordinates: y, x Data variables: pwat (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 cloudBase float64 0.0 latitude (y, x) float64 12.19 12.39 12.58 12.77 ... 57.68 57.49 57.29 longitude (y, x) float64 226.5 227.2 227.9 228.7 ... 308.5 309.6 310.6 valid_time datetime64[ns] 2018-09-17 Dimensions without coordinates: y, x Data variables: pres (y, x) float32 ... gh (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 cloudTop float64 0.0 latitude (y, x) float64 12.19 12.39 12.58 12.77 ... 57.68 57.49 57.29 longitude (y, x) float64 226.5 227.2 227.9 228.7 ... 308.5 309.6 310.6 valid_time datetime64[ns] 2018-09-17 Dimensions without coordinates: y, x Data variables: pres (y, x) float32 ... t (y, x) float32 ... gh (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 heightAboveGround float64 10.0 latitude (y, x) float64 ... longitude (y, x) float64 ... valid_time datetime64[ns] ... Dimensions without coordinates: y, x Data variables: u10 (y, x) float32 ... v10 (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 heightAboveGround float64 2.0 latitude (y, x) float64 12.19 12.39 12.58 ... 57.68 57.49 57.29 longitude (y, x) float64 226.5 227.2 227.9 ... 308.5 309.6 310.6 valid_time datetime64[ns] 2018-09-17 Dimensions without coordinates: y, x Data variables: t2m (y, x) float32 ... r2 (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (heightAboveGroundLayer: 2, y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 * heightAboveGroundLayer (heightAboveGroundLayer) float64 1e+03 3e+03 latitude (y, x) float64 ... longitude (y, x) float64 ... valid_time datetime64[ns] ... Dimensions without coordinates: y, x Data variables: hlcy (heightAboveGroundLayer, y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (isobaricInhPa: 19, y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 * isobaricInhPa (isobaricInhPa) float64 1e+03 950.0 900.0 ... 150.0 100.0 latitude (y, x) float64 12.19 12.39 12.58 12.77 ... 57.68 57.49 57.29 longitude (y, x) float64 226.5 227.2 227.9 228.7 ... 308.5 309.6 310.6 valid_time datetime64[ns] 2018-09-17 Dimensions without coordinates: y, x Data variables: t (isobaricInhPa, y, x) float32 ... u (isobaricInhPa, y, x) float32 ... v (isobaricInhPa, y, x) float32 ... w (isobaricInhPa, y, x) float32 ... gh (isobaricInhPa, y, x) float32 ... r (isobaricInhPa, y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (isobaricInhPa: 5, y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 * isobaricInhPa (isobaricInhPa) float64 1e+03 850.0 700.0 500.0 250.0 latitude (y, x) float64 ... longitude (y, x) float64 ... valid_time datetime64[ns] ... Dimensions without coordinates: y, x Data variables: absv (isobaricInhPa, y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 isothermZero float64 0.0 latitude (y, x) float64 12.19 12.39 12.58 12.77 ... 57.68 57.49 57.29 longitude (y, x) float64 226.5 227.2 227.9 228.7 ... 308.5 309.6 310.6 valid_time datetime64[ns] 2018-09-17 Dimensions without coordinates: y, x Data variables: gh (y, x) float32 ... r (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 maxWind float64 0.0 latitude (y, x) float64 12.19 12.39 12.58 12.77 ... 57.68 57.49 57.29 longitude (y, x) float64 226.5 227.2 227.9 228.7 ... 308.5 309.6 310.6 valid_time datetime64[ns] 2018-09-17 Dimensions without coordinates: y, x Data variables: pres (y, x) float32 ... u (y, x) float32 ... v (y, x) float32 ... gh (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 meanSea float64 0.0 latitude (y, x) float64 12.19 12.39 12.58 12.77 ... 57.68 57.49 57.29 longitude (y, x) float64 226.5 227.2 227.9 228.7 ... 308.5 309.6 310.6 valid_time datetime64[ns] 2018-09-17 Dimensions without coordinates: y, x Data variables: prmsl (y, x) float32 ... mslet (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (pressureFromGroundLayer: 2, y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 * pressureFromGroundLayer (pressureFromGroundLayer) float64 9e+03 1.8e+04 latitude (y, x) float64 12.19 12.39 12.58 ... 57.49 57.29 longitude (y, x) float64 226.5 227.2 227.9 ... 309.6 310.6 valid_time datetime64[ns] 2018-09-17 Dimensions without coordinates: y, x Data variables: cape (pressureFromGroundLayer, y, x) float32 ... cin (pressureFromGroundLayer, y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (pressureFromGroundLayer: 5, y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 * pressureFromGroundLayer (pressureFromGroundLayer) float64 3e+03 ... 1.5e+04 latitude (y, x) float64 12.19 12.39 12.58 ... 57.49 57.29 longitude (y, x) float64 226.5 227.2 227.9 ... 309.6 310.6 valid_time datetime64[ns] 2018-09-17 Dimensions without coordinates: y, x Data variables: t (pressureFromGroundLayer, y, x) float32 ... u (pressureFromGroundLayer, y, x) float32 ... v (pressureFromGroundLayer, y, x) float32 ... r (pressureFromGroundLayer, y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 pressureFromGroundLayer float64 3e+03 latitude (y, x) float64 ... longitude (y, x) float64 ... valid_time datetime64[ns] ... Dimensions without coordinates: y, x Data variables: pli (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 pressureFromGroundLayer float64 1.8e+04 latitude (y, x) float64 ... longitude (y, x) float64 ... valid_time datetime64[ns] ... Dimensions without coordinates: y, x Data variables: 4lftx (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 surface float64 0.0 latitude (y, x) float64 12.19 12.39 12.58 12.77 ... 57.68 57.49 57.29 longitude (y, x) float64 226.5 227.2 227.9 228.7 ... 308.5 309.6 310.6 valid_time datetime64[ns] 2018-09-17 Dimensions without coordinates: y, x Data variables: cape (y, x) float32 ... sp (y, x) float32 ... acpcp (y, x) float32 ... cin (y, x) float32 ... orog (y, x) float32 ... tp (y, x) float32 ... crain (y, x) float32 ... cfrzr (y, x) float32 ... cicep (y, x) float32 ... csnow (y, x) float32 ... gust (y, x) float32 ... hpbl (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP , <xarray.Dataset> Dimensions: (y: 65, x: 93) Coordinates: time datetime64[ns] 2018-09-17 step timedelta64[ns] 00:00:00 tropopause float64 0.0 latitude (y, x) float64 12.19 12.39 12.58 12.77 ... 57.68 57.49 57.29 longitude (y, x) float64 226.5 227.2 227.9 228.7 ... 308.5 309.6 310.6 valid_time datetime64[ns] 2018-09-17 Dimensions without coordinates: y, x Data variables: t (y, x) float32 ... u (y, x) float32 ... v (y, x) float32 ... trpp (y, x) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: kwbc GRIB_centreDescription: US National Weather Service - NCEP... GRIB_subCentre: 0 Conventions: CF-1.7 institution: US National Weather Service - NCEP ]
Advanced usage
Write support
Please note that write support is Alpha. Only xarray.Dataset’s in canonical form, that is, with the coordinates names matching exactly the cfgrib coordinates, can be saved at the moment:
>>> from cfgrib.xarray_to_grib import to_grib >>> ds = xr.open_dataset('era5-levels-members.grib', engine='cfgrib').sel(number=0) >>> ds <xarray.Dataset> Dimensions: (time: 4, isobaricInhPa: 2, latitude: 61, longitude: 120) Coordinates: number int64 0 * time (time) datetime64[ns] 2017-01-01 ... 2017-01-02T12:00:00 step timedelta64[ns] ... * isobaricInhPa (isobaricInhPa) float64 850.0 500.0 * latitude (latitude) float64 90.0 87.0 84.0 81.0 ... -84.0 -87.0 -90.0 * longitude (longitude) float64 0.0 3.0 6.0 9.0 ... 351.0 354.0 357.0 valid_time (time) datetime64[ns] ... Data variables: z (time, isobaricInhPa, latitude, longitude) float32 ... t (time, isobaricInhPa, latitude, longitude) float32 ... Attributes: GRIB_edition: 1 GRIB_centre: ecmf GRIB_centreDescription: European Centre for Medium-Range Weather Forecasts GRIB_subCentre: 0 Conventions: CF-1.7 institution: European Centre for Medium-Range Weather Forecasts history: ... >>> to_grib(ds, 'out1.grib', grib_keys={'edition': 2}) >>> xr.open_dataset('out1.grib', engine='cfgrib') <xarray.Dataset> Dimensions: (time: 4, isobaricInhPa: 2, latitude: 61, longitude: 120) Coordinates: number ... * time (time) datetime64[ns] 2017-01-01 ... 2017-01-02T12:00:00 step timedelta64[ns] ... * isobaricInhPa (isobaricInhPa) float64 850.0 500.0 * latitude (latitude) float64 90.0 87.0 84.0 81.0 ... -84.0 -87.0 -90.0 * longitude (longitude) float64 0.0 3.0 6.0 9.0 ... 351.0 354.0 357.0 valid_time (time) datetime64[ns] ... Data variables: z (time, isobaricInhPa, latitude, longitude) float32 ... t (time, isobaricInhPa, latitude, longitude) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: ecmf GRIB_centreDescription: European Centre for Medium-Range Weather Forecasts GRIB_subCentre: 0 Conventions: CF-1.7 institution: European Centre for Medium-Range Weather Forecasts history: ...
Per-variable GRIB keys can be set by setting the attrs variable with key prefixed by GRIB_, for example:
>>> import numpy as np >>> import xarray as xr >>> ds2 = xr.DataArray( ... np.zeros((5, 6)) + 300., ... coords=[ ... np.linspace(90., -90., 5), ... np.linspace(0., 360., 6, endpoint=False), ... ], ... dims=['latitude', 'longitude'], ... ).to_dataset(name='skin_temperature') >>> ds2.skin_temperature.attrs['GRIB_shortName'] = 'skt' >>> to_grib(ds2, 'out2.grib') >>> xr.open_dataset('out2.grib', engine='cfgrib') <xarray.Dataset> Dimensions: (latitude: 5, longitude: 6) Coordinates: time datetime64[ns] ... step timedelta64[ns] ... surface float64 ... * latitude (latitude) float64 90.0 45.0 0.0 -45.0 -90.0 * longitude (longitude) float64 0.0 60.0 120.0 180.0 240.0 300.0 valid_time datetime64[ns] ... Data variables: skt (latitude, longitude) float32 ... Attributes: GRIB_edition: 2 GRIB_centre: consensus GRIB_centreDescription: Consensus GRIB_subCentre: 0 Conventions: CF-1.7 institution: Consensus history: ...
Dataset / Variable API
The use of xarray is not mandatory and you can access the content of a GRIB file as an hypercube with the high level API in a Python interpreter:
>>> ds = cfgrib.open_file('era5-levels-members.grib') >>> ds.attributes['GRIB_edition'] 1 >>> sorted(ds.dimensions.items()) [('isobaricInhPa', 2), ('latitude', 61), ('longitude', 120), ('number', 10), ('time', 4)] >>> sorted(ds.variables) ['isobaricInhPa', 'latitude', 'longitude', 'number', 'step', 't', 'time', 'valid_time', 'z'] >>> var = ds.variables['t'] >>> var.dimensions ('number', 'time', 'isobaricInhPa', 'latitude', 'longitude') >>> var.data[:, :, :, :, :].mean() 262.92133 >>> ds = cfgrib.open_file('era5-levels-members.grib') >>> ds.attributes['GRIB_edition'] 1 >>> sorted(ds.dimensions.items()) [('isobaricInhPa', 2), ('latitude', 61), ('longitude', 120), ('number', 10), ('time', 4)] >>> sorted(ds.variables) ['isobaricInhPa', 'latitude', 'longitude', 'number', 'step', 't', 'time', 'valid_time', 'z'] >>> var = ds.variables['t'] >>> var.dimensions ('number', 'time', 'isobaricInhPa', 'latitude', 'longitude') >>> var.data[:, :, :, :, :].mean() 262.92133
GRIB index file
By default cfgrib saves the index of the GRIB file to disk appending .idx to the GRIB file name. Index files are an experimental and completely optional feature, feel free to remove them and try again in case of problems. Index files saving can be disable passing adding indexpath='' to the backend_kwargs keyword argument.
Project resources
Development |
|
Download |
|
User support |
|
Code quality |
Contributing
The main repository is hosted on GitHub, testing, bug reports and contributions are highly welcomed and appreciated:
https://github.com/ecmwf/cfgrib
Please see the CONTRIBUTING.rst document for the best way to help.
Lead developers:
Baudouin Raoult - ECMWF
Main contributors:
Aureliana Barghini - B-Open
Leonardo Barcaroli - B-Open
See also the list of contributors who participated in this project.
License
Copyright 2017-2021 European Centre for Medium-Range Weather Forecasts (ECMWF).
Licensed under the Apache License, Version 2.0 (the “License”); you may not use this file except in compliance with the License. You may obtain a copy of the License at: http://www.apache.org/licenses/LICENSE-2.0. Unless required by applicable law or agreed to in writing, software distributed under the License is distributed on an “AS IS” BASIS, WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. See the License for the specific language governing permissions and limitations under the License.
Changelog for cfgrib
0.9.10.2 (unreleased)
Nothing changed yet.
0.9.10.1 (2022-03-16)
Fix failure to read index files. See #292.
Allow backend kwargs to be provided in the to_netcdf executable, either via a json format string, or a path to a json file via -b. See #288.
Fixed issue where the use of relpath() could cause a problem on Windows. See #284.
Fix passing of pathlib.Path. See #282.
Fixed issue where writing an ensemble number into a GRIB file caused an error. See #278.
0.9.10.0 (2022-01-31)
Big internal refactor to add support for a generic Fieldset similar to Metview. See #243.
0.9.9.1 (2021-09-29)
0.9.9.0 (2021-04-09)
Depend on the ECMWF eccodes python package to access the low level ecCodes C-library, dropping all other GRIB decoding options. See: #95, #14. #204, #147 and #141.
Many performance improvements during the generation of the index and during data access. See: #142 and #197.
filter_by_keys now can select on all keys known to ecCodes without the need to add non default ones to read_keys explicitly. See: #187.
Include support for engine=”cfgrib” using xarray 0.18+ new backend API. See: #216.
Fixed issue where could not load a GRIB message that has only one grid point. See: #199.
Decode level coordinates as float in all cases, fixed issue with non-int levels. See: #195.
0.9.8.5 (2020-11-11)
Simpler and clearer messages in the event of errors.
Use ECCODES_DIR environment variable if present. Ported from eccodes-python by xavierabellan. See: #162.
Fix using current ecCodes bindings when setting CFGRIB_USE_EXTERNAL_ECCODES_BINDINGS=1.
0.9.8.4 (2020-08-03)
Use ecmwflibs if present to find the ecCodes installation.
0.9.8.3 (2020-06-25)
Added support for indexingDate, indexingTime time coordinates.
lambert_azimuthal_equal_area grids are now returned as 2D arrays. See: #119.
0.9.8.2 (2020-05-22)
0.9.8.1 (2020-03-13)
Always open GRIB files in binary mode, by @b8raoult
0.9.8.0 (2020-03-12)
Add support of experimental pyeccodes low-level driver by @b8raoult
0.9.7.7 (2020-01-24)
Add support for forecastMonth in cf2cdm.translate_coords.
0.9.7.6 (2019-12-05)
Fix the README.
0.9.7.5 (2019-12-05)
Deprecate ensure_valid_time and the config option preferred_time_dimension that are now better handled via time_dims.
0.9.7.4 (2019-11-22)
Add more options to time_dims forecasts products may be represented as ('time', 'verifying_time') or ('time', 'forecastMonth'). See: #97.
0.9.7.3 (2019-11-04)
Add support for selecting the time coordinates to use as dimensions via time_dims. Forecasts products may be represented as ('time', 'step') (the default), ('time', 'valid_time') or ('valid_time', 'step'). See: #97.
Reduce the in-memory footprint of the FieldIndex and the size of .idx files.
0.9.7.2 (2019-09-24)
0.9.7.1 (2019-07-08)
0.9.7 (2019-05-27)
0.9.7rc1 (2019-05-14)
Drop support for Python 2, in line with xarray 0.12.0. The 0.9.6.x series will be supported long term for Python 2 users. See: #69.
Sync internal ecCodes bindings API to the one in eccodes-python. See: #81.
Source code has been formatted with black -S -l 99.
Added initial support for spectral coordinates.
0.9.6.2 (2019-04-15)
Improve merging of variables into a dataset. See: #63.
0.9.6.1.post1 (2019-03-17)
Fix an issue in the README format.
0.9.6.1 (2019-03-17)
Fixed (for real) MULTI-FIELD messages, See: #45.
Added a protocol version to the index file. Old *.idx files must be removed.
0.9.6.post1 (2019-03-07)
Fix an important typo in the README. See: #64.
0.9.6 (2019-02-26)
0.9.5.7 (2019-02-24)
Fixed a serious bug in the computation of the suggested filter_by_keys for non-cubic GRIB files. As a result cfgrib.xarray_store.open_datasets was not finding all the variables in the files. See: #54.
Fixed a serious bug in variable naming that could drop or at worse mix the values of variables. Again see: #54.
Re-opened #45 as the fix was returning wrong data. Now we are back to dropping all variable in a MULTI-FIELD except the first.
0.9.5.6 (2019-02-04)
Do not set explicit timezone in units to avoid crashing some versions of xarray. See: #44.
0.9.5.5 (2019-02-02)
Enable ecCodes implicit MULTI-FIELD support by default, needed for NAM Products by NCEP. See: #45.
Added support for depthBelowLand coordinate.
0.9.5.4 (2019-01-25)
Add support for building valid_time from a bad time-step hypercube.
0.9.5.3 (2019-01-25)
Also convert is valid_time can index all times and steps in translate_coords.
0.9.5.2 (2019-01-24)
Set valid_time as preferred time dimension for the CDS data model.
Fall back to using the generic GRIB2 ecCodes template when no better option is found. See: #39.
0.9.5.1 (2018-12-27)
0.9.5 (2018-12-20)
Drop support for xarray versions prior to v0.11 to reduce complexity. (This is really only v0.10.9). See: #32.
Declare the data as CF-1.7 compliant via the Conventions global attribute. See: #36.
Tested larger-than-memory and distributed processing via dask and dask.distributed. See: #33.
Promote write support via cfgrib.to_grib to Alpha. See: #18.
Provide the cf2cdm.translate_coords utility function to translate the coordinates between CF-compliant data models, defined by out_name, units and store_direction. See: #24.
Provide cfgrib.__version__. See: #31.
Raise with a better error message when users attempt to open a file that is not a GRIB. See: #34.
Make 2D grids for rotated_ll and rotated_gg gridType’s. See: #35.
0.9.4.1 (2018-11-08)
Fix formatting for PyPI page.
0.9.4 (2018-11-08)
Saves one index file per set of index_keys in a much more robust way.
Refactor CF-encoding and add the new encode_cf option to backend_kwargs. See: #23.
Refactor error handling and the option to ignore errors (not well documented yet). See: #13.
Do not crash on gridType not fully supported by the installed ecCodes See: #27.
Several smaller bug fixes and performance improvements.
0.9.3.1 (2018-10-28)
Assorted README fixes, in particular advertise index file support as alpha.
0.9.3 (2018-10-28)
Big performance improvement: add alpha support to save to and read from disk the GRIB index produced by the full-file scan at the first open. See: #20.
0.9.2 (2018-10-22)
Rename coordinate air_pressure to isobaricInhPa for consistency with all other vertical level coordinates. See: #25.
0.9.1.post1 (2018-10-19)
Fix PyPI description.
0.9.1 (2018-10-19)
Change the usage of cfgrib.open_dataset to allign it with xarray.open_dataset, in particular filter_by_key must be added into the backend_kwargs dictionary. See: #21.
0.9.0 (2018-10-14)
Beta release with read support.
Project details
Release history Release notifications | RSS feed
Download files
Download the file for your platform. If you're not sure which to choose, learn more about installing packages.
Source Distribution
Built Distribution
File details
Details for the file cfgrib-0.9.10.1.tar.gz
.
File metadata
- Download URL: cfgrib-0.9.10.1.tar.gz
- Upload date:
- Size: 6.4 MB
- Tags: Source
- Uploaded using Trusted Publishing? No
- Uploaded via: twine/3.1.1 pkginfo/1.5.0.1 requests/2.25.1 setuptools/56.0.0 requests-toolbelt/0.9.1 tqdm/4.40.0 CPython/3.6.8
File hashes
Algorithm | Hash digest | |
---|---|---|
SHA256 | 34ef4a0334c1077f6ad7c15de1d8a79a0aa0980187b2848deb8cb8acb3755b74 |
|
MD5 | 1c3650aacd8cfbd9c5d35cfc03307b31 |
|
BLAKE2b-256 | 8195673e98a5ca4bc70060cc042e59baad78f95a466ed7e2d3ac0143cf1f9a4a |
File details
Details for the file cfgrib-0.9.10.1-py3-none-any.whl
.
File metadata
- Download URL: cfgrib-0.9.10.1-py3-none-any.whl
- Upload date:
- Size: 46.0 kB
- Tags: Python 3
- Uploaded using Trusted Publishing? No
- Uploaded via: twine/3.1.1 pkginfo/1.5.0.1 requests/2.25.1 setuptools/56.0.0 requests-toolbelt/0.9.1 tqdm/4.40.0 CPython/3.6.8
File hashes
Algorithm | Hash digest | |
---|---|---|
SHA256 | f16615de9e8b81b55b267a7bcbc22db470065c7e7989c79c7d56c47b22d9bd40 |
|
MD5 | 1ea6a9892efed23d3b652a728bb1f6fa |
|
BLAKE2b-256 | b0c50eb15b1a648c43104957ef73ae72d8ee27cc94b67cf935ca0e3d31416209 |