grib2io

Introduction

grib2io is a Python package that provides an interface to the NCEP GRIB2 C (g2c) library for the purpose of reading and writing WMO GRIdded Binary, Edition 2 (GRIB2) messages. A physical file can contain one or more GRIB2 messages.

GRIB2 file IO is performed directly in Python. The unpacking/packing of GRIB2 integer, coded metadata and data sections is performed by the g2c library functions via the g2clib Cython extension module. The decoding/encoding of GRIB2 metadata is translated into more descriptive, plain language metadata by providing a mapping of the integer coded metadata to the appropriate GRIB2 code tables. The NCEP GRIB2 code tables are a part of the grib2io package.

Index File

As of v2.6.0, grib2io provides the capability to create (and read) an index file for a GRIB2 file. The index file contents are specific to Python and grib2io in that pickle is used to dump the index dictionary to file. Using index files can dramatically improve performance in situations where the same file will be read multiple times. The index file name is the original GRIB2 file name with a hash string appended, followed by the file extension, .grib2ioidx. The hash string is the SHA-1 of the GRIB2 file name and the file size. For example, GRIB2 file,

gfs.t00z.pgrb2.1p00.f024

when opened, grib2io will generate an index file with the following name,

gfs.t00z.pgrb2.1p00.f024_0422a93bfd6d095bd0a942ba5e9fe42e76050123.grib2ioidx

By default, grib2io will generate an index file or use an existing one. The generation and usage of a grib2io index file can be turned off by providing kwargs use_index=False and/or save_index=False in grib2io.open().

Interpolation

As of v2.4.0, spatial interpolation via NCEPLIPS-ip Fortran library is now a part of the grib2io package. Previously, interpolation was handled by a separate component package, grib2io-interp, which is now deprecated. grib2io-interp provided interpolation via a F2PY-generated interface to NCEPLIBS-ip, which has become difficult since the removal of distutils from Python 3.12+.

NCEPLIBS-ip interpolation Fortran subroutines contain the BIND(C) attribute which provides an equivalent C-interface. grib2io now provides a Cython-based interface, iplib, to these Fortran subroutines via their C-interface. If NCEPLIBS-ip was built with OpenMP support, iplib will provide functions for getting and setting the number of OpenMP threads.

Xarray Backend

grib2io provides a Xarray backend engine so that many GRIB2 messages can be represented as N-dimensional DataArray objects and collected along common coordinates as Datasets or DataTrees. The Xarray backend engine API is experimental and is subject to change without backward compatibility.

Tutorials

The following Jupyter Notebooks are available as tutorials:

 1from ._grib2io import (
 2    open,
 3    interpolate,
 4    interpolate_to_stations,
 5    Grib2Message,
 6    _Grib2Message,
 7    Grib2GridDef,
 8    msgs_from_index,
 9    __doc__,
10)
11from . import tables, templates, utils
12
13try:
14    from . import __config__
15
16    __version__ = __config__.grib2io_version
17    has_interpolation = __config__.has_interpolation
18    has_openmp_support = __config__.has_openmp_support
19    g2c_static = __config__.g2c_static
20    ip_static = __config__.ip_static
21    extra_objects = __config__.extra_objects
22except ImportError:
23    pass
24
25from .g2clib import __version__ as __g2clib_version__
26from .g2clib import _has_jpeg
27from .g2clib import _has_png
28from .g2clib import _has_aec
29
30from .tables.originating_centers import _ncep_grib2_table_version
31
32__all__ = [
33    "open",
34    "show_config",
35    "interpolate",
36    "interpolate_to_stations",
37    "tables",
38    "templates",
39    "utils",
40    "Grib2Message",
41    "_Grib2Message",
42    "Grib2GridDef",
43    "msgs_from_index",
44    "__doc__",
45]
46
47has_jpeg_support = bool(_has_jpeg)
48has_png_support = bool(_has_png)
49has_aec_support = bool(_has_aec)
50
51ncep_grib2_table_version = _ncep_grib2_table_version
52g2c_version = __g2clib_version__
53
54
55def show_config():
56    """Print grib2io build configuration information."""
57    print(f"grib2io version {__version__} Configuration:")
58    print("")
59    print(f"NCEPLIBS-g2c library version: {__g2clib_version__}")
60    print(f"\tStatic library: {g2c_static}")
61    print(f"\tJPEG compression support: {has_jpeg_support}")
62    print(f"\tPNG compression support: {has_png_support}")
63    print(f"\tAEC compression support: {has_aec_support}")
64    print("")
65    print(f"NCEPLIPS-ip support: {has_interpolation}")
66    print(f"\tStatic library: {ip_static}")
67    print(f"\tOpenMP support: {has_openmp_support}")
68    print("")
69    print("Static libs:")
70    for lib in extra_objects:
71        print(f"\t{lib}")
72    print("")
73    print(f"NCEP GRIB2 Table Version: {_ncep_grib2_table_version}")
class open:
117class open:
118    """
119    GRIB2 File Object.
120
121    This class can accommodate a physical file with one or more GRIB2
122    messages or a bytes object containing a GRIB2 messages.
123
124    A physical file can contain one or more GRIB2 messages.  When instantiated,
125    class `grib2io.open`, the file named `filename` is opened for reading (`mode
126    = 'r'`) and is automatically indexed.  The indexing procedure reads some of
127    the GRIB2 metadata for all GRIB2 Messages.  A GRIB2 Message may contain
128    submessages whereby Section 2-7 can be repeated.  grib2io accommodates for
129    this by flattening any GRIB2 submessages into multiple individual messages.
130
131    It is important to note that GRIB2 files from some Meteorological agencies
132    contain other data than GRIB2 messages.  GRIB2 files from ECMWF can contain
133    GRIB1 and GRIB2 messages.  grib2io checks for these and safely ignores them.
134
135    Attributes
136    ----------
137    closed : bool
138        `True` is file handle is close; `False` otherwise.
139    current_message : int
140        Current position of the file in units of GRIB2 Messages. (read only)
141    indexfile : str
142        Index file for the GRIB2 file.
143    levels : tuple
144        Tuple containing a unique list of wgrib2-formatted level/layer strings.
145    messages : int
146        Count of GRIB2 Messages contained in the file.
147    mode : str
148        File IO mode of opening the file.
149    name : str
150        Full path name of the GRIB2 file.
151    save_index : bool
152        Whether to save a pickle-based index file for the GRIB2 file. Default is `True`.
153    size : int
154        Size of the file in units of bytes.
155    use_index
156        Whether to use an existing pickle-based index file for the GRIB2 file. Default is `True`.
157    variables : tuple
158        Tuple containing a unique list of variable short names (i.e. GRIB2
159        abbreviation names).
160    """
161
162    __slots__ = (
163        "_fileid",
164        "_filehandle",
165        "_hasindex",
166        "_index",
167        "_msgs",
168        "_pos",
169        "closed",
170        "current_message",
171        "indexfile",
172        "messages",
173        "mode",
174        "name",
175        "size",
176        "save_index",
177        "use_index",
178    )
179
180    def __init__(
181        self,
182        filename: Union[bytes, str, Path, IO],
183        mode: Literal["r", "w", "x"] = "r",
184        *,
185        save_index=True,
186        use_index=True,
187        _xarray_backend=False,
188        **kwargs,
189    ):
190        """
191        Initialize GRIB2 File object instance.
192
193        Parameters
194        ----------
195        filename
196            File path containing GRIB2 messages OR bytes OR file-like object.
197        mode
198            File access mode where "r" opens the files for reading only; "w"
199            opens the file for overwriting and "x" for writing to a new file.
200        save_index
201            Whether to save a pickle-based index file for the GRIB2 file. Default is True.
202
203            .. versionadded:: 2.6.0
204        use_index
205            Whether to use an existing pickle-based index file for the GRIB2 file. Default is True.
206
207            .. versionadded:: 2.6.0
208        """
209
210        # All write modes are read/write; all modes are binary.
211        if mode in ("a", "x", "w"):
212            mode += "+"
213        mode = mode + "b"
214
215        self._hasindex = False
216        self.indexfile = None
217        self.mode = mode
218        self.save_index = save_index
219        self.size = 0
220        self.use_index = use_index
221
222        if isinstance(filename, bytes):
223            if "r" not in self.mode:
224                raise ValueError(
225                    "Invalid mode for bytes input: GRIB2 data supplied as bytes is read-only. Use mode='r' or provide a filename instead."
226                )
227
228            self.current_message = 0
229            if filename[:2] == _GZIP_HEADER:
230                filename = gzip.decompress(filename)
231            if filename[:4] + filename[-4:] != b"GRIB7777":
232                raise ValueError("Invalid GRIB bytes")
233            self._filehandle = BytesIO(filename)
234            self.name = "<in-memory-file>"
235            self.size = len(filename)
236            self._fileid = hashlib.sha1(
237                (self.name + str(self.size)).encode("ASCII")
238            ).hexdigest()
239            self._index = build_index(self._filehandle)
240            self._msgs = msgs_from_index(self._index, filehandle=self._filehandle)
241            self.messages = len(self._msgs)
242
243        elif all(hasattr(filename, attr) for attr in ("read", "seek", "tell")):
244            # Handle file-like objects (including S3File, etc.)
245            self.current_message = 0
246            self._filehandle = filename
247            self.name = getattr(filename, "name", filename.__repr__())
248            # For file-like objects, we can't get file stats, so create a unique ID
249            self._fileid = hashlib.sha1(
250                (self.name + str(id(filename))).encode("ASCII")
251            ).hexdigest()
252            # Disable index file usage for file-like objects since we can't save/load them
253            self.use_index = False
254            self.save_index = False
255
256            if "r" in self.mode:
257                # Build index from the file-like object
258                self._index = build_index(self._filehandle)
259                self._msgs = msgs_from_index(self._index, filehandle=self._filehandle)
260                self.messages = len(self._msgs)
261
262            # Get size for s3fs.core.S3File object.
263            self.size = self._filehandle.info()["size"] or 0
264
265        else:
266            self.current_message = 0
267            if "r" in mode:
268                self._filehandle = builtins.open(filename, mode=mode)
269                # Some GRIB2 files are gzipped, so check for that here, but
270                # raise error when using xarray backend.
271                # Gzip files contain a 2-byte header b'\x1f\x8b'.
272                if self._filehandle.read(2) == _GZIP_HEADER:
273                    self._filehandle.close()
274                    if _xarray_backend:
275                        raise RuntimeError(
276                            "Gzip GRIB2 files are not supported by the Xarray backend."
277                        )
278                    self._filehandle = gzip.open(filename, mode=mode)
279                else:
280                    self._filehandle = builtins.open(filename, mode=mode)
281            else:
282                self._filehandle = builtins.open(filename, mode=mode)
283                self.name = os.path.abspath(filename)
284            self.name = os.path.abspath(filename)
285            fstat = os.stat(self.name)
286            self.size = fstat.st_size
287            self._fileid = hashlib.sha1(
288                (self.name + str(fstat.st_ino) + str(self.size)).encode("ASCII")
289            ).hexdigest()
290
291            if "r" in self.mode:
292                self.indexfile = f"{self.name}_{self._fileid}.grib2ioidx"
293                if self.use_index and os.path.exists(self.indexfile):
294                    try:
295                        with builtins.open(self.indexfile, "rb") as file:
296                            index = pickle.load(file)
297                        self._index = index
298                    except Exception as e:
299                        warnings.warn(
300                            f"found indexfile: {self.indexfile}, but unable to load it: {e}\n"
301                            f"re-forming index from grib2file, but not writing indexfile"
302                        )
303                        self._index = build_index(self._filehandle)
304                    self._hasindex = True
305                else:
306                    self._index = build_index(self._filehandle)
307                    if self.save_index:
308                        try:
309                            serialize_index(self._index, self.indexfile)
310                        except Exception as e:
311                            warnings.warn(
312                                f"index was not serialized for future use: {e}"
313                            )
314
315                self._msgs = msgs_from_index(self._index, filehandle=self._filehandle)
316
317                self.messages = len(self._msgs)
318            elif "w" or "x" in self.mode:
319                self.messages = 0
320                self.current_message = None
321
322        self.closed = self._filehandle.closed
323
324    def __delete__(self, instance):
325        self.close()
326        del self._index
327
328    def __enter__(self):
329        return self
330
331    def __exit__(self, atype, value, traceback):
332        self.close()
333
334    def __iter__(self):
335        yield from self._msgs
336
337    def __len__(self):
338        return self.messages
339
340    def __repr__(self):
341        strings = []
342        for k in self.__slots__:
343            if k.startswith("_"):
344                continue
345            strings.append("%s = %s\n" % (k, eval("self." + k)))
346        return "".join(strings)
347
348    def __getitem__(self, key):
349        if isinstance(key, int):
350            if abs(key) >= len(self._msgs):
351                raise IndexError("index out of range")
352            else:
353                return self._msgs[key]
354        elif isinstance(key, str):
355            return self.select(shortName=key)
356        elif isinstance(key, slice):
357            return self._msgs[key]
358        elif isinstance(key, (list, tuple, set)):
359            if len(key) == 0:
360                return iter(self._msgs)
361            indices = sorted(key) if isinstance(key, set) else key
362
363            def _iter_msgs():
364                for i in indices:
365                    if not isinstance(i, int):
366                        raise TypeError(
367                            f"indices must be integers; got {type(i).__name__}"
368                        )
369                    if abs(i) >= len(self._msgs):
370                        raise IndexError(f"index out of range: {i}")
371                    yield self._msgs[i]
372
373            return _iter_msgs()
374        else:
375            raise KeyError(
376                "Key must be an integer, slice, GRIB2 variable shortName, "
377                "or an iterable of integer indices."
378            )
379
380    @property
381    def levels(self):
382        """Provides a unique tuple of level strings."""
383        if self._hasindex:
384            return tuple(sorted(set([msg.level for msg in self._msgs])))
385        else:
386            return None
387
388    @property
389    def variables(self):
390        """Provides a unique tuple of variable shortName strings."""
391        if self._hasindex:
392            return tuple(sorted(set([msg.shortName for msg in self._msgs])))
393        else:
394            return None
395
396    def close(self):
397        """Close the file handle."""
398        if not self._filehandle.closed:
399            self.messages = 0
400            self.current_message = 0
401            self._filehandle.close()
402            self.closed = self._filehandle.closed
403
404    def read(self, size: Optional[int] = None):
405        """
406        Read size amount of GRIB2 messages from the current position.
407
408        If no argument is given, then size is None and all messages are returned
409        from the current position in the file. This read method follows the
410        behavior of Python's builtin open() function, but whereas that operates
411        on units of bytes, we operate on units of GRIB2 messages.
412
413        Parameters
414        ----------
415        size: default=None
416            The number of GRIB2 messages to read from the current position. If
417            no argument is give, the default value is None and remainder of
418            the file is read.
419
420        Returns
421        -------
422        read
423            ``Grib2Message`` object when size = 1 or a list of Grib2Messages
424            when size > 1.
425        """
426        if size is not None and size < 0:
427            size = None
428        if size is None or size > 1:
429            start = self.tell()
430            stop = self.messages if size is None else start + size
431            if size is None:
432                self.current_message = self.messages - 1
433            else:
434                self.current_message += size
435            return self._msgs[slice(start, stop, 1)]
436        elif size == 1:
437            self.current_message += 1
438            return self._msgs[self.current_message]
439        else:
440            None
441
442    def seek(self, pos: int):
443        """
444        Set the position within the file in units of GRIB2 messages.
445
446        Parameters
447        ----------
448        pos
449            The GRIB2 Message number to set the file pointer to.
450        """
451        if self._hasindex:
452            self._filehandle.seek(self._index["sectionOffset"][0][pos])
453            self.current_message = pos
454
455    def tell(self):
456        """Returns the position of the file in units of GRIB2 Messages."""
457        return self.current_message
458
459    def select(self, **kwargs):
460        """Select GRIB2 messages by `Grib2Message` attributes."""
461        # TODO: Added ability to process multiple values for each keyword (attribute)
462        idxs = []
463        nkeys = len(kwargs.keys())
464        for k, v in kwargs.items():
465            for m in self._msgs:
466                if hasattr(m, k) and getattr(m, k) == v:
467                    idxs.append(m._msgnum)
468        idxs = np.array(idxs, dtype=np.int32)
469        return [
470            self._msgs[i]
471            for i in [
472                ii[0]
473                for ii in collections.Counter(idxs).most_common()
474                if ii[1] == nkeys
475            ]
476        ]
477
478    def write(self, msg):
479        """
480        Writes GRIB2 message object to file.
481
482        Parameters
483        ----------
484        msg
485            GRIB2 message objects to write to file.
486        """
487        if isinstance(msg, list):
488            for m in msg:
489                self.write(m)
490            return
491
492        if issubclass(msg.__class__, _Grib2Message):
493            # TODO: We can consider letting pack return packed bytes instead of associating with message object
494            if hasattr(msg, "_msg"):
495                # write already packed bytes
496                self._filehandle.write(msg._msg)
497            else:
498                if (
499                    msg._signature == msg._generate_signature()
500                    and msg._data is None
501                    and hasattr(msg._ondiskarray, "filehandle")
502                ):
503                    # write unchanged message from input
504                    offset = msg._ondiskarray.filehandle.tell()
505                    msg._ondiskarray.filehandle.seek(msg._ondiskarray.offset)
506                    self._filehandle.write(
507                        msg._ondiskarray.filehandle.read(msg.section0[-1])
508                    )
509                    msg._ondiskarray.filehandle.seek(offset)
510                else:
511                    msg.pack()
512                    self._filehandle.write(msg._msg)
513            self.size = os.path.getsize(self.name)
514            self.messages += 1
515        else:
516            raise TypeError("msg must be a Grib2Message object.")
517        return
518
519    def flush(self):
520        """Flush the file object buffer."""
521        self._filehandle.flush()
522
523    def levels_by_var(self, name: str):
524        """
525        Return a list of level strings given a variable shortName.
526
527        Parameters
528        ----------
529        name
530            Grib2Message variable shortName
531
532        Returns
533        -------
534        levels_by_var
535            A list of unique level strings.
536        """
537        return list(sorted(set([msg.level for msg in self.select(shortName=name)])))
538
539    def vars_by_level(self, level: str):
540        """
541        Return a list of variable shortName strings given a level.
542
543        Parameters
544        ----------
545        level
546            Grib2Message variable level
547
548        Returns
549        -------
550        vars_by_level
551            A list of unique variable shortName strings.
552        """
553        return list(sorted(set([msg.shortName for msg in self.select(level=level)])))

GRIB2 File Object.

This class can accommodate a physical file with one or more GRIB2 messages or a bytes object containing a GRIB2 messages.

A physical file can contain one or more GRIB2 messages. When instantiated, class grib2io.open, the file named filename is opened for reading (mode = 'r') and is automatically indexed. The indexing procedure reads some of the GRIB2 metadata for all GRIB2 Messages. A GRIB2 Message may contain submessages whereby Section 2-7 can be repeated. grib2io accommodates for this by flattening any GRIB2 submessages into multiple individual messages.

It is important to note that GRIB2 files from some Meteorological agencies contain other data than GRIB2 messages. GRIB2 files from ECMWF can contain GRIB1 and GRIB2 messages. grib2io checks for these and safely ignores them.

Attributes
  • closed (bool): True is file handle is close; False otherwise.
  • current_message (int): Current position of the file in units of GRIB2 Messages. (read only)
  • indexfile (str): Index file for the GRIB2 file.
  • levels (tuple): Tuple containing a unique list of wgrib2-formatted level/layer strings.
  • messages (int): Count of GRIB2 Messages contained in the file.
  • mode (str): File IO mode of opening the file.
  • name (str): Full path name of the GRIB2 file.
  • save_index (bool): Whether to save a pickle-based index file for the GRIB2 file. Default is True.
  • size (int): Size of the file in units of bytes.
  • use_index: Whether to use an existing pickle-based index file for the GRIB2 file. Default is True.
  • variables (tuple): Tuple containing a unique list of variable short names (i.e. GRIB2 abbreviation names).
open( filename: Union[bytes, str, pathlib.Path, IO], mode: Literal['r', 'w', 'x'] = 'r', *, save_index=True, use_index=True, _xarray_backend=False, **kwargs)
180    def __init__(
181        self,
182        filename: Union[bytes, str, Path, IO],
183        mode: Literal["r", "w", "x"] = "r",
184        *,
185        save_index=True,
186        use_index=True,
187        _xarray_backend=False,
188        **kwargs,
189    ):
190        """
191        Initialize GRIB2 File object instance.
192
193        Parameters
194        ----------
195        filename
196            File path containing GRIB2 messages OR bytes OR file-like object.
197        mode
198            File access mode where "r" opens the files for reading only; "w"
199            opens the file for overwriting and "x" for writing to a new file.
200        save_index
201            Whether to save a pickle-based index file for the GRIB2 file. Default is True.
202
203            .. versionadded:: 2.6.0
204        use_index
205            Whether to use an existing pickle-based index file for the GRIB2 file. Default is True.
206
207            .. versionadded:: 2.6.0
208        """
209
210        # All write modes are read/write; all modes are binary.
211        if mode in ("a", "x", "w"):
212            mode += "+"
213        mode = mode + "b"
214
215        self._hasindex = False
216        self.indexfile = None
217        self.mode = mode
218        self.save_index = save_index
219        self.size = 0
220        self.use_index = use_index
221
222        if isinstance(filename, bytes):
223            if "r" not in self.mode:
224                raise ValueError(
225                    "Invalid mode for bytes input: GRIB2 data supplied as bytes is read-only. Use mode='r' or provide a filename instead."
226                )
227
228            self.current_message = 0
229            if filename[:2] == _GZIP_HEADER:
230                filename = gzip.decompress(filename)
231            if filename[:4] + filename[-4:] != b"GRIB7777":
232                raise ValueError("Invalid GRIB bytes")
233            self._filehandle = BytesIO(filename)
234            self.name = "<in-memory-file>"
235            self.size = len(filename)
236            self._fileid = hashlib.sha1(
237                (self.name + str(self.size)).encode("ASCII")
238            ).hexdigest()
239            self._index = build_index(self._filehandle)
240            self._msgs = msgs_from_index(self._index, filehandle=self._filehandle)
241            self.messages = len(self._msgs)
242
243        elif all(hasattr(filename, attr) for attr in ("read", "seek", "tell")):
244            # Handle file-like objects (including S3File, etc.)
245            self.current_message = 0
246            self._filehandle = filename
247            self.name = getattr(filename, "name", filename.__repr__())
248            # For file-like objects, we can't get file stats, so create a unique ID
249            self._fileid = hashlib.sha1(
250                (self.name + str(id(filename))).encode("ASCII")
251            ).hexdigest()
252            # Disable index file usage for file-like objects since we can't save/load them
253            self.use_index = False
254            self.save_index = False
255
256            if "r" in self.mode:
257                # Build index from the file-like object
258                self._index = build_index(self._filehandle)
259                self._msgs = msgs_from_index(self._index, filehandle=self._filehandle)
260                self.messages = len(self._msgs)
261
262            # Get size for s3fs.core.S3File object.
263            self.size = self._filehandle.info()["size"] or 0
264
265        else:
266            self.current_message = 0
267            if "r" in mode:
268                self._filehandle = builtins.open(filename, mode=mode)
269                # Some GRIB2 files are gzipped, so check for that here, but
270                # raise error when using xarray backend.
271                # Gzip files contain a 2-byte header b'\x1f\x8b'.
272                if self._filehandle.read(2) == _GZIP_HEADER:
273                    self._filehandle.close()
274                    if _xarray_backend:
275                        raise RuntimeError(
276                            "Gzip GRIB2 files are not supported by the Xarray backend."
277                        )
278                    self._filehandle = gzip.open(filename, mode=mode)
279                else:
280                    self._filehandle = builtins.open(filename, mode=mode)
281            else:
282                self._filehandle = builtins.open(filename, mode=mode)
283                self.name = os.path.abspath(filename)
284            self.name = os.path.abspath(filename)
285            fstat = os.stat(self.name)
286            self.size = fstat.st_size
287            self._fileid = hashlib.sha1(
288                (self.name + str(fstat.st_ino) + str(self.size)).encode("ASCII")
289            ).hexdigest()
290
291            if "r" in self.mode:
292                self.indexfile = f"{self.name}_{self._fileid}.grib2ioidx"
293                if self.use_index and os.path.exists(self.indexfile):
294                    try:
295                        with builtins.open(self.indexfile, "rb") as file:
296                            index = pickle.load(file)
297                        self._index = index
298                    except Exception as e:
299                        warnings.warn(
300                            f"found indexfile: {self.indexfile}, but unable to load it: {e}\n"
301                            f"re-forming index from grib2file, but not writing indexfile"
302                        )
303                        self._index = build_index(self._filehandle)
304                    self._hasindex = True
305                else:
306                    self._index = build_index(self._filehandle)
307                    if self.save_index:
308                        try:
309                            serialize_index(self._index, self.indexfile)
310                        except Exception as e:
311                            warnings.warn(
312                                f"index was not serialized for future use: {e}"
313                            )
314
315                self._msgs = msgs_from_index(self._index, filehandle=self._filehandle)
316
317                self.messages = len(self._msgs)
318            elif "w" or "x" in self.mode:
319                self.messages = 0
320                self.current_message = None
321
322        self.closed = self._filehandle.closed

Initialize GRIB2 File object instance.

Parameters
  • filename: File path containing GRIB2 messages OR bytes OR file-like object.
  • mode: File access mode where "r" opens the files for reading only; "w" opens the file for overwriting and "x" for writing to a new file.
  • save_index: Whether to save a pickle-based index file for the GRIB2 file. Default is True.

New in version 2.6.0.

  • use_index: Whether to use an existing pickle-based index file for the GRIB2 file. Default is True.

New in version 2.6.0.

indexfile
mode
save_index
size
use_index
closed
levels
380    @property
381    def levels(self):
382        """Provides a unique tuple of level strings."""
383        if self._hasindex:
384            return tuple(sorted(set([msg.level for msg in self._msgs])))
385        else:
386            return None

Provides a unique tuple of level strings.

variables
388    @property
389    def variables(self):
390        """Provides a unique tuple of variable shortName strings."""
391        if self._hasindex:
392            return tuple(sorted(set([msg.shortName for msg in self._msgs])))
393        else:
394            return None

Provides a unique tuple of variable shortName strings.

def close(self):
396    def close(self):
397        """Close the file handle."""
398        if not self._filehandle.closed:
399            self.messages = 0
400            self.current_message = 0
401            self._filehandle.close()
402            self.closed = self._filehandle.closed

Close the file handle.

def read(self, size: Optional[int] = None):
404    def read(self, size: Optional[int] = None):
405        """
406        Read size amount of GRIB2 messages from the current position.
407
408        If no argument is given, then size is None and all messages are returned
409        from the current position in the file. This read method follows the
410        behavior of Python's builtin open() function, but whereas that operates
411        on units of bytes, we operate on units of GRIB2 messages.
412
413        Parameters
414        ----------
415        size: default=None
416            The number of GRIB2 messages to read from the current position. If
417            no argument is give, the default value is None and remainder of
418            the file is read.
419
420        Returns
421        -------
422        read
423            ``Grib2Message`` object when size = 1 or a list of Grib2Messages
424            when size > 1.
425        """
426        if size is not None and size < 0:
427            size = None
428        if size is None or size > 1:
429            start = self.tell()
430            stop = self.messages if size is None else start + size
431            if size is None:
432                self.current_message = self.messages - 1
433            else:
434                self.current_message += size
435            return self._msgs[slice(start, stop, 1)]
436        elif size == 1:
437            self.current_message += 1
438            return self._msgs[self.current_message]
439        else:
440            None

Read size amount of GRIB2 messages from the current position.

If no argument is given, then size is None and all messages are returned from the current position in the file. This read method follows the behavior of Python's builtin open() function, but whereas that operates on units of bytes, we operate on units of GRIB2 messages.

Parameters
  • size (default=None): The number of GRIB2 messages to read from the current position. If no argument is give, the default value is None and remainder of the file is read.
Returns
  • read: Grib2Message object when size = 1 or a list of Grib2Messages when size > 1.
def seek(self, pos: int):
442    def seek(self, pos: int):
443        """
444        Set the position within the file in units of GRIB2 messages.
445
446        Parameters
447        ----------
448        pos
449            The GRIB2 Message number to set the file pointer to.
450        """
451        if self._hasindex:
452            self._filehandle.seek(self._index["sectionOffset"][0][pos])
453            self.current_message = pos

Set the position within the file in units of GRIB2 messages.

Parameters
  • pos: The GRIB2 Message number to set the file pointer to.
def tell(self):
455    def tell(self):
456        """Returns the position of the file in units of GRIB2 Messages."""
457        return self.current_message

Returns the position of the file in units of GRIB2 Messages.

def select(self, **kwargs):
459    def select(self, **kwargs):
460        """Select GRIB2 messages by `Grib2Message` attributes."""
461        # TODO: Added ability to process multiple values for each keyword (attribute)
462        idxs = []
463        nkeys = len(kwargs.keys())
464        for k, v in kwargs.items():
465            for m in self._msgs:
466                if hasattr(m, k) and getattr(m, k) == v:
467                    idxs.append(m._msgnum)
468        idxs = np.array(idxs, dtype=np.int32)
469        return [
470            self._msgs[i]
471            for i in [
472                ii[0]
473                for ii in collections.Counter(idxs).most_common()
474                if ii[1] == nkeys
475            ]
476        ]

Select GRIB2 messages by Grib2Message attributes.

def write(self, msg):
478    def write(self, msg):
479        """
480        Writes GRIB2 message object to file.
481
482        Parameters
483        ----------
484        msg
485            GRIB2 message objects to write to file.
486        """
487        if isinstance(msg, list):
488            for m in msg:
489                self.write(m)
490            return
491
492        if issubclass(msg.__class__, _Grib2Message):
493            # TODO: We can consider letting pack return packed bytes instead of associating with message object
494            if hasattr(msg, "_msg"):
495                # write already packed bytes
496                self._filehandle.write(msg._msg)
497            else:
498                if (
499                    msg._signature == msg._generate_signature()
500                    and msg._data is None
501                    and hasattr(msg._ondiskarray, "filehandle")
502                ):
503                    # write unchanged message from input
504                    offset = msg._ondiskarray.filehandle.tell()
505                    msg._ondiskarray.filehandle.seek(msg._ondiskarray.offset)
506                    self._filehandle.write(
507                        msg._ondiskarray.filehandle.read(msg.section0[-1])
508                    )
509                    msg._ondiskarray.filehandle.seek(offset)
510                else:
511                    msg.pack()
512                    self._filehandle.write(msg._msg)
513            self.size = os.path.getsize(self.name)
514            self.messages += 1
515        else:
516            raise TypeError("msg must be a Grib2Message object.")
517        return

Writes GRIB2 message object to file.

Parameters
  • msg: GRIB2 message objects to write to file.
def flush(self):
519    def flush(self):
520        """Flush the file object buffer."""
521        self._filehandle.flush()

Flush the file object buffer.

def levels_by_var(self, name: str):
523    def levels_by_var(self, name: str):
524        """
525        Return a list of level strings given a variable shortName.
526
527        Parameters
528        ----------
529        name
530            Grib2Message variable shortName
531
532        Returns
533        -------
534        levels_by_var
535            A list of unique level strings.
536        """
537        return list(sorted(set([msg.level for msg in self.select(shortName=name)])))

Return a list of level strings given a variable shortName.

Parameters
  • name: Grib2Message variable shortName
Returns
  • levels_by_var: A list of unique level strings.
def vars_by_level(self, level: str):
539    def vars_by_level(self, level: str):
540        """
541        Return a list of variable shortName strings given a level.
542
543        Parameters
544        ----------
545        level
546            Grib2Message variable level
547
548        Returns
549        -------
550        vars_by_level
551            A list of unique variable shortName strings.
552        """
553        return list(sorted(set([msg.shortName for msg in self.select(level=level)])))

Return a list of variable shortName strings given a level.

Parameters
  • level: Grib2Message variable level
Returns
  • vars_by_level: A list of unique variable shortName strings.
current_message
messages
name
def show_config():
56def show_config():
57    """Print grib2io build configuration information."""
58    print(f"grib2io version {__version__} Configuration:")
59    print("")
60    print(f"NCEPLIBS-g2c library version: {__g2clib_version__}")
61    print(f"\tStatic library: {g2c_static}")
62    print(f"\tJPEG compression support: {has_jpeg_support}")
63    print(f"\tPNG compression support: {has_png_support}")
64    print(f"\tAEC compression support: {has_aec_support}")
65    print("")
66    print(f"NCEPLIPS-ip support: {has_interpolation}")
67    print(f"\tStatic library: {ip_static}")
68    print(f"\tOpenMP support: {has_openmp_support}")
69    print("")
70    print("Static libs:")
71    for lib in extra_objects:
72        print(f"\t{lib}")
73    print("")
74    print(f"NCEP GRIB2 Table Version: {_ncep_grib2_table_version}")

Print grib2io build configuration information.

def interpolate( a, method: Union[int, str], grid_def_in, grid_def_out, method_options=None, num_threads=1):
2231def interpolate(
2232    a,
2233    method: Union[int, str],
2234    grid_def_in,
2235    grid_def_out,
2236    method_options=None,
2237    num_threads=1,
2238):
2239    """
2240    This is the module-level interpolation function.
2241
2242    This interfaces with the 4-byte (32-bit float) [NCEPLIBS-ip library](https://github.com/NOAA-EMC/NCEPLIBS-ip)
2243    through grib2io's internal iplib Cython extension module.
2244
2245    Parameters
2246    ----------
2247    a : numpy.ndarray or tuple
2248        Input data.  If `a` is a `numpy.ndarray`, scalar interpolation will be
2249        performed.  If `a` is a `tuple`, then vector interpolation will be
2250        performed with the assumption that u = a[0] and v = a[1] and are both
2251        `numpy.ndarray`.
2252
2253        These data are expected to be in 2-dimensional form with shape (ny, nx)
2254        or 3-dimensional (:, ny, nx) where the 1st dimension represents another
2255        spatial, temporal, or classification (i.e. ensemble members) dimension.
2256        The function will properly flatten the (ny,nx) dimensions into (nx * ny)
2257        acceptable for input into the interpolation subroutines. If needed, these
2258        data will be converted to `np.float32`.
2259    method
2260        Interpolate method to use. This can either be an integer or string using
2261        the following mapping:
2262
2263        | Interpolate Scheme | Integer Value |
2264        | :---:              | :---:         |
2265        | 'bilinear'         | 0             |
2266        | 'bicubic'          | 1             |
2267        | 'neighbor'         | 2             |
2268        | 'budget'           | 3             |
2269        | 'spectral'         | 4             |
2270        | 'neighbor-budget'  | 6             |
2271
2272    grid_def_in : grib2io.Grib2GridDef
2273        Grib2GridDef object for the input grid.
2274    grid_def_out : grib2io.Grib2GridDef
2275        Grib2GridDef object for the output grid or station points.
2276    method_options : list of ints, optional
2277        Interpolation options. See the NCEPLIBS-ip documentation for
2278        more information on how these are used.
2279    num_threads : int, optional
2280        Number of OpenMP threads to use for interpolation. The default
2281        value is 1. If NCEPLIBS-ip and grib2io's iplib extension module
2282        was not built with OpenMP, then this keyword argument and value
2283        will have no impact.
2284
2285    Returns
2286    -------
2287    interpolate
2288        Returns a `numpy.ndarray` of dtype `np.float32` when scalar interpolation
2289        is performed or a `tuple` of `numpy.ndarray`s when vector interpolation is
2290        performed with the assumptions that 0-index is the interpolated u and
2291        1-index is the interpolated v.
2292    """
2293
2294    try:
2295        from . import iplib
2296    except ImportError:
2297        raise ImportError(
2298            "NCEPLIBS-ip library not found. Interpolation is not available."
2299        )
2300
2301    prev_num_threads = 1
2302    try:
2303        prev_num_threads = iplib.openmp_get_num_threads()
2304        iplib.openmp_set_num_threads(num_threads)
2305    except AttributeError:
2306        pass
2307
2308    print(
2309        f"grib2io.interpolate thread report: OpenMP num threads = {iplib.openmp_get_num_threads()}"
2310    )
2311
2312    if isinstance(method, int) and method not in _interp_schemes.values():
2313        raise ValueError("Invalid interpolation method.")
2314    elif isinstance(method, str):
2315        if method in _interp_schemes.keys():
2316            method = _interp_schemes[method]
2317        else:
2318            raise ValueError("Invalid interpolation method.")
2319
2320    if method_options is None:
2321        method_options = np.zeros((20), dtype=np.int32)
2322        if method in {3, 6}:
2323            method_options[0:2] = -1
2324
2325    mi = grid_def_in.npoints
2326    mo = grid_def_out.npoints
2327
2328    # Adjust shape of input array(s)
2329    a, newshp = _adjust_array_shape_for_interp(a, grid_def_in, grid_def_out)
2330
2331    # Call interpolation subroutines according to type of a.
2332    if isinstance(a, np.ndarray):
2333        # Scalar
2334        km = a.shape[0]
2335        if np.any(np.isnan(a)):
2336            ibi = np.ones((km), dtype=np.int32)
2337            li = np.where(np.isnan(a), 0, 1).astype(np.uint8)
2338        else:
2339            ibi = np.zeros((km), dtype=np.int32)
2340            li = np.zeros(a.shape, dtype=np.uint8)
2341        no, rlat, rlon, ibo, lo, go, iret = iplib.interpolate_scalar(
2342            method,
2343            method_options,
2344            grid_def_in.gdtn,
2345            np.array(grid_def_in.gdt, dtype=np.int32),
2346            grid_def_out.gdtn,
2347            np.array(grid_def_out.gdt, dtype=np.int32),
2348            mi,
2349            mo,
2350            km,
2351            ibi,
2352            li,
2353            a,
2354        )
2355        out = np.where(lo == 0, np.nan, go).reshape(newshp)
2356    elif isinstance(a, tuple):
2357        # Vector
2358        km = a[0].shape[0]
2359        if np.any(np.isnan(a)):
2360            ibi = np.ones((km), dtype=np.int32)
2361            li = np.where(np.isnan(a), 0, 1).astype(np.uint8)
2362        else:
2363            ibi = np.zeros((km), dtype=np.int32)
2364            li = np.zeros(a[0].shape, dtype=np.uint8)
2365        no, rlat, rlon, crot, srot, ibo, lo, uo, vo, iret = iplib.interpolate_vector(
2366            method,
2367            method_options,
2368            grid_def_in.gdtn,
2369            np.array(grid_def_in.gdt, dtype=np.int32),
2370            grid_def_out.gdtn,
2371            np.array(grid_def_out.gdt, dtype=np.int32),
2372            mi,
2373            mo,
2374            km,
2375            ibi,
2376            li,
2377            a[0],
2378            a[1],
2379        )
2380        uo = np.where(lo == 0, np.nan, uo).reshape(newshp)
2381        vo = np.where(lo == 0, np.nan, vo).reshape(newshp)
2382        out = (uo, vo)
2383
2384    try:
2385        iplib.openmp_set_num_threads(prev_num_threads)
2386    except AttributeError:
2387        pass
2388
2389    return out

This is the module-level interpolation function.

This interfaces with the 4-byte (32-bit float) NCEPLIBS-ip library through grib2io's internal iplib Cython extension module.

Parameters
  • a (numpy.ndarray or tuple): Input data. If a is a numpy.ndarray, scalar interpolation will be performed. If a is a tuple, then vector interpolation will be performed with the assumption that u = a[0] and v = a[1] and are both numpy.ndarray.

    These data are expected to be in 2-dimensional form with shape (ny, nx) or 3-dimensional (:, ny, nx) where the 1st dimension represents another spatial, temporal, or classification (i.e. ensemble members) dimension. The function will properly flatten the (ny,nx) dimensions into (nx * ny) acceptable for input into the interpolation subroutines. If needed, these data will be converted to np.float32.

  • method: Interpolate method to use. This can either be an integer or string using the following mapping:
Interpolate Scheme Integer Value
'bilinear' 0
'bicubic' 1
'neighbor' 2
'budget' 3
'spectral' 4
'neighbor-budget' 6

  • grid_def_in (grib2io.Grib2GridDef): Grib2GridDef object for the input grid.
  • grid_def_out (grib2io.Grib2GridDef): Grib2GridDef object for the output grid or station points.
  • method_options (list of ints, optional): Interpolation options. See the NCEPLIBS-ip documentation for more information on how these are used.
  • num_threads (int, optional): Number of OpenMP threads to use for interpolation. The default value is 1. If NCEPLIBS-ip and grib2io's iplib extension module was not built with OpenMP, then this keyword argument and value will have no impact.
Returns
  • interpolate: Returns a numpy.ndarray of dtype np.float32 when scalar interpolation is performed or a tuple of numpy.ndarrays when vector interpolation is performed with the assumptions that 0-index is the interpolated u and 1-index is the interpolated v.
def interpolate_to_stations( a, method: Union[int, str], grid_def_in, lats, lons, method_options=None, num_threads=1):
2392def interpolate_to_stations(
2393    a,
2394    method: Union[int, str],
2395    grid_def_in,
2396    lats,
2397    lons,
2398    method_options=None,
2399    num_threads=1,
2400):
2401    """
2402    Module-level interpolation function for interpolation to stations.
2403
2404    Interfaces with the 4-byte (32-bit float) [NCEPLIBS-ip library](https://github.com/NOAA-EMC/NCEPLIBS-ip)
2405    via grib2io's iplib Cython exntension module. It supports scalar and
2406    vector interpolation according to the type of object a.
2407
2408    Parameters
2409    ----------
2410    a : numpy.ndarray or tuple
2411        Input data.  If `a` is a `numpy.ndarray`, scalar interpolation will be
2412        performed.  If `a` is a `tuple`, then vector interpolation will be
2413        performed with the assumption that u = a[0] and v = a[1] and are both
2414        `numpy.ndarray`.
2415
2416        These data are expected to be in 2-dimensional form with shape (ny, nx)
2417        or 3-dimensional (:, ny, nx) where the 1st dimension represents another
2418        spatial, temporal, or classification (i.e. ensemble members) dimension.
2419        The function will properly flatten the (ny,nx) dimensions into (nx * ny)
2420        acceptable for input into the interpolation subroutines. If needed, these
2421        data will be converted to `np.float32`.
2422    method
2423        Interpolate method to use. This can either be an integer or string using
2424        the following mapping:
2425
2426        | Interpolate Scheme | Integer Value |
2427        | :---:              | :---:         |
2428        | 'bilinear'         | 0             |
2429        | 'bicubic'          | 1             |
2430        | 'neighbor'         | 2             |
2431        | 'budget'           | 3             |
2432        | 'spectral'         | 4             |
2433        | 'neighbor-budget'  | 6             |
2434
2435    grid_def_in : grib2io.Grib2GridDef
2436        Grib2GridDef object for the input grid.
2437    lats : numpy.ndarray or list
2438        Latitudes for station points
2439    lons : numpy.ndarray or list
2440        Longitudes for station points
2441    method_options : list of ints, optional
2442        Interpolation options. See the NCEPLIBS-ip documentation for
2443        more information on how these are used.
2444    num_threads : int, optional
2445        Number of OpenMP threads to use for interpolation. The default
2446        value is 1. If NCEPLIBS-ip and grib2io's iplib extension module
2447        was not built with OpenMP, then this keyword argument and value
2448        will have no impact.
2449
2450    Returns
2451    -------
2452    interpolate_to_stations
2453        Returns a `numpy.ndarray` of dtype `np.float32` when scalar
2454        interpolation is performed or a `tuple` of `numpy.ndarray`s
2455        when vector interpolation is performed with the assumptions
2456        that 0-index is the interpolated u and 1-index is the
2457        interpolated v.
2458    """
2459    try:
2460        from . import iplib
2461    except ImportError:
2462        raise ImportError(
2463            "NCEPLIBS-ip library not found. Interpolation is not available."
2464        )
2465
2466    # Define function to apply mask when stations are outside grid domain
2467    def _reshape_and_mask_post_interp(a, shape, mask):
2468        a = a.reshape(shape)
2469        if a.shape[-1] != mask.shape[0]:
2470            raise ValueError("Station mask length does not match interpolated data.")
2471        a[..., mask] = np.nan
2472        return a
2473
2474    prev_num_threads = 1
2475    try:
2476        prev_num_threads = iplib.openmp_get_num_threads()
2477        iplib.openmp_set_num_threads(num_threads)
2478    except AttributeError:
2479        pass
2480
2481    if isinstance(method, int) and method not in _interp_schemes.values():
2482        raise ValueError("Invalid interpolation method.")
2483    elif isinstance(method, str):
2484        if method in _interp_schemes.keys():
2485            method = _interp_schemes[method]
2486        else:
2487            raise ValueError("Invalid interpolation method.")
2488
2489    if method_options is None:
2490        method_options = np.zeros((20), dtype=np.int32)
2491        if method in {3, 6}:
2492            method_options[0:2] = -1
2493
2494    # Check lats and lons
2495    if isinstance(lats, list):
2496        nlats = len(lats)
2497    elif isinstance(lats, np.ndarray) and len(lats.shape) == 1:
2498        nlats = lats.shape[0]
2499    else:
2500        raise ValueError("Station latitudes must be a list or 1-D NumPy array.")
2501    if isinstance(lons, list):
2502        nlons = len(lons)
2503    elif isinstance(lons, np.ndarray) and len(lons.shape) == 1:
2504        nlons = lons.shape[0]
2505    else:
2506        raise ValueError("Station longitudes must be a list or 1-D NumPy array.")
2507    if nlats != nlons:
2508        raise ValueError("Station lats and lons must be same size.")
2509
2510    mi = grid_def_in.npoints
2511    mo = nlats
2512
2513    # Adjust shape of input array(s)
2514    a, newshp = _adjust_array_shape_for_interp_stations(a, grid_def_in, mo)
2515
2516    # Use gdtn = -1 for stations and an empty template array
2517    gdtn_out = -1
2518    gdt_out = np.zeros((200), dtype=np.int32)
2519
2520    # Before we interpolate, get the grid coordinates for stations.
2521    xloc, yloc = utils.latlon_to_ij(
2522        grid_def_in.gdtn,
2523        grid_def_in.gdt,
2524        np.array(lats, dtype=np.float32),
2525        np.array(lons, dtype=np.float32),
2526    )
2527    xloc_mask = np.where(np.isnan(xloc), True, False)
2528    yloc_mask = np.where(np.isnan(yloc), True, False)
2529    mask = xloc_mask & yloc_mask
2530
2531    # Call interpolation subroutines according to type of a.
2532    if isinstance(a, np.ndarray):
2533        # Scalar
2534        km = a.shape[0]
2535        ibi = np.zeros((km), dtype=np.int32)
2536        li = np.zeros(a.shape, dtype=np.uint8)
2537        no, rlat, rlon, ibo, lo, go, iret = iplib.interpolate_scalar(
2538            method,
2539            method_options,
2540            grid_def_in.gdtn,
2541            np.array(grid_def_in.gdt, dtype=np.int32),
2542            gdtn_out,
2543            gdt_out,
2544            mi,
2545            mo,
2546            km,
2547            ibi,
2548            li,
2549            a,
2550            lats=np.array(lats, dtype=np.float32),
2551            lons=np.array(lons, dtype=np.float32),
2552        )
2553        out = _reshape_and_mask_post_interp(go, newshp, mask)
2554
2555    elif isinstance(a, tuple):
2556        # Vector
2557        km = a[0].shape[0]
2558        ibi = np.zeros((km), dtype=np.int32)
2559        li = np.zeros(a[0].shape, dtype=np.uint8)
2560        no, rlat, rlon, crot, srot, ibo, lo, uo, vo, iret = iplib.interpolate_vector(
2561            method,
2562            method_options,
2563            grid_def_in.gdtn,
2564            np.array(grid_def_in.gdt, dtype=np.int32),
2565            gdtn_out,
2566            gdt_out,
2567            mi,
2568            mo,
2569            km,
2570            ibi,
2571            li,
2572            a[0],
2573            a[1],
2574            lats=np.array(lats, dtype=np.float32),
2575            lons=np.array(lons, dtype=np.float32),
2576        )
2577        out = (
2578            _reshape_and_mask_post_interp(uo, newshp, mask),
2579            _reshape_and_mask_post_interp(vo, newshp, mask),
2580        )
2581
2582    try:
2583        iplib.openmp_set_num_threads(prev_num_threads)
2584    except AttributeError:
2585        pass
2586
2587    return out

Module-level interpolation function for interpolation to stations.

Interfaces with the 4-byte (32-bit float) NCEPLIBS-ip library via grib2io's iplib Cython exntension module. It supports scalar and vector interpolation according to the type of object a.

Parameters
  • a (numpy.ndarray or tuple): Input data. If a is a numpy.ndarray, scalar interpolation will be performed. If a is a tuple, then vector interpolation will be performed with the assumption that u = a[0] and v = a[1] and are both numpy.ndarray.

    These data are expected to be in 2-dimensional form with shape (ny, nx) or 3-dimensional (:, ny, nx) where the 1st dimension represents another spatial, temporal, or classification (i.e. ensemble members) dimension. The function will properly flatten the (ny,nx) dimensions into (nx * ny) acceptable for input into the interpolation subroutines. If needed, these data will be converted to np.float32.

  • method: Interpolate method to use. This can either be an integer or string using the following mapping:
Interpolate Scheme Integer Value
'bilinear' 0
'bicubic' 1
'neighbor' 2
'budget' 3
'spectral' 4
'neighbor-budget' 6

  • grid_def_in (grib2io.Grib2GridDef): Grib2GridDef object for the input grid.
  • lats (numpy.ndarray or list): Latitudes for station points
  • lons (numpy.ndarray or list): Longitudes for station points
  • method_options (list of ints, optional): Interpolation options. See the NCEPLIBS-ip documentation for more information on how these are used.
  • num_threads (int, optional): Number of OpenMP threads to use for interpolation. The default value is 1. If NCEPLIBS-ip and grib2io's iplib extension module was not built with OpenMP, then this keyword argument and value will have no impact.
Returns
  • interpolate_to_stations: Returns a numpy.ndarray of dtype np.float32 when scalar interpolation is performed or a tuple of numpy.ndarrays when vector interpolation is performed with the assumptions that 0-index is the interpolated u and 1-index is the interpolated v.
class Grib2Message:
 903class Grib2Message:
 904    """
 905    Creation class for a GRIB2 message.
 906
 907    This class returns a dynamically-created Grib2Message object that
 908    inherits from `_Grib2Message` and grid, product, data representation
 909    template classes according to the template numbers for the respective
 910    sections. If `section3`, `section4`, or `section5` are omitted, then
 911    the appropriate keyword arguments for the template number `gdtn=`,
 912    `pdtn=`, or `drtn=` must be provided.
 913
 914    Parameters
 915    ----------
 916    section0
 917        GRIB2 section 0 array.
 918    section1
 919        GRIB2 section 1 array.
 920    section2
 921        Local Use section data.
 922    section3
 923        GRIB2 section 3 array.
 924    section4
 925        GRIB2 section 4 array.
 926    section5
 927        GRIB2 section 5 array.
 928
 929    Returns
 930    -------
 931    Msg
 932        A dynamically-create Grib2Message object that inherits from
 933        _Grib2Message, a grid definition template class, product
 934        definition template class, and a data representation template
 935        class.
 936    """
 937
 938    def __new__(
 939        self,
 940        section0: NDArray = np.array([struct.unpack(">I", b"GRIB")[0], 0, 0, 2, 0]),
 941        section1: NDArray = np.zeros((13), dtype=np.int64),
 942        section2: Optional[bytes] = None,
 943        section3: Optional[NDArray] = None,
 944        section4: Optional[NDArray] = None,
 945        section5: Optional[NDArray] = None,
 946        *args,
 947        **kwargs,
 948    ):
 949
 950        if np.all(section1 == 0):
 951            try:
 952                # Python >= 3.10
 953                section1[5:11] = datetime.datetime.fromtimestamp(
 954                    0, datetime.UTC
 955                ).timetuple()[:6]
 956            except AttributeError:
 957                # Python < 3.10
 958                section1[5:11] = datetime.datetime.utcfromtimestamp(0).timetuple()[:6]
 959
 960        bases = list()
 961        if section3 is None:
 962            if "gdtn" in kwargs.keys():
 963                gdtn = kwargs["gdtn"]
 964                Gdt = templates.gdt_class_by_gdtn(gdtn)
 965                bases.append(Gdt)
 966                section3 = np.zeros((Gdt._len + 5), dtype=np.int64)
 967                section3[4] = gdtn
 968            else:
 969                raise ValueError(
 970                    "Must provide GRIB2 Grid Definition Template Number or section 3 array"
 971                )
 972        else:
 973            gdtn = section3[4]
 974            Gdt = templates.gdt_class_by_gdtn(gdtn)
 975            bases.append(Gdt)
 976
 977        if section4 is None:
 978            if "pdtn" in kwargs.keys():
 979                pdtn = kwargs["pdtn"]
 980                Pdt = templates.pdt_class_by_pdtn(pdtn)
 981                bases.append(Pdt)
 982                section4 = np.zeros((Pdt._len + 2), dtype=np.int64)
 983                section4[1] = pdtn
 984            else:
 985                raise ValueError(
 986                    "Must provide GRIB2 Production Definition Template Number or section 4 array"
 987                )
 988        else:
 989            pdtn = section4[1]
 990            Pdt = templates.pdt_class_by_pdtn(pdtn)
 991            bases.append(Pdt)
 992
 993        if section5 is None:
 994            if "drtn" in kwargs.keys():
 995                drtn = kwargs["drtn"]
 996                Drt = templates.drt_class_by_drtn(drtn)
 997                bases.append(Drt)
 998                section5 = np.zeros((Drt._len + 2), dtype=np.int64)
 999                section5[1] = drtn
1000            else:
1001                raise ValueError(
1002                    "Must provide GRIB2 Data Representation Template Number or section 5 array"
1003                )
1004        else:
1005            drtn = section5[1]
1006            Drt = templates.drt_class_by_drtn(drtn)
1007            bases.append(Drt)
1008
1009        # attempt to use existing Msg class if it has already been made with gdtn,pdtn,drtn combo
1010        try:
1011            Msg = _msg_class_store[f"{gdtn}:{pdtn}:{drtn}"]
1012        except KeyError:
1013
1014            @dataclass(init=False, repr=False)
1015            class Msg(_Grib2Message, *bases):
1016                pass
1017
1018            _msg_class_store[f"{gdtn}:{pdtn}:{drtn}"] = Msg
1019
1020        return Msg(section0, section1, section2, section3, section4, section5, *args)

Creation class for a GRIB2 message.

This class returns a dynamically-created Grib2Message object that inherits from _Grib2Message and grid, product, data representation template classes according to the template numbers for the respective sections. If section3, section4, or section5 are omitted, then the appropriate keyword arguments for the template number gdtn=, pdtn=, or drtn= must be provided.

Parameters
  • section0: GRIB2 section 0 array.
  • section1: GRIB2 section 1 array.
  • section2: Local Use section data.
  • section3: GRIB2 section 3 array.
  • section4: GRIB2 section 4 array.
  • section5: GRIB2 section 5 array.
Returns
  • Msg: A dynamically-create Grib2Message object that inherits from _Grib2Message, a grid definition template class, product definition template class, and a data representation template class.
@dataclass
class _Grib2Message:
1023@dataclass
1024class _Grib2Message:
1025    """
1026    GRIB2 Message base class.
1027    """
1028
1029    # GRIB2 Sections
1030    section0: NDArray = field(init=True, repr=False)
1031    section1: NDArray = field(init=True, repr=False)
1032    section2: bytes = field(init=True, repr=False)
1033    section3: NDArray = field(init=True, repr=False)
1034    section4: NDArray = field(init=True, repr=False)
1035    section5: NDArray = field(init=True, repr=False)
1036    bitMapFlag: templates.Grib2Metadata = field(init=True, repr=False, default=255)
1037
1038    # Section 0 looked up attributes
1039    indicatorSection: NDArray = field(
1040        init=False, repr=False, default=templates.IndicatorSection()
1041    )
1042    discipline: templates.Grib2Metadata = field(
1043        init=False, repr=False, default=templates.Discipline()
1044    )
1045
1046    # Section 1 looked up attributes
1047    identificationSection: NDArray = field(
1048        init=False, repr=False, default=templates.IdentificationSection()
1049    )
1050    originatingCenter: templates.Grib2Metadata = field(
1051        init=False, repr=False, default=templates.OriginatingCenter()
1052    )
1053    originatingSubCenter: templates.Grib2Metadata = field(
1054        init=False, repr=False, default=templates.OriginatingSubCenter()
1055    )
1056    masterTableInfo: templates.Grib2Metadata = field(
1057        init=False, repr=False, default=templates.MasterTableInfo()
1058    )
1059    localTableInfo: templates.Grib2Metadata = field(
1060        init=False, repr=False, default=templates.LocalTableInfo()
1061    )
1062    significanceOfReferenceTime: templates.Grib2Metadata = field(
1063        init=False, repr=False, default=templates.SignificanceOfReferenceTime()
1064    )
1065    year: int = field(init=False, repr=False, default=templates.Year())
1066    month: int = field(init=False, repr=False, default=templates.Month())
1067    day: int = field(init=False, repr=False, default=templates.Day())
1068    hour: int = field(init=False, repr=False, default=templates.Hour())
1069    minute: int = field(init=False, repr=False, default=templates.Minute())
1070    second: int = field(init=False, repr=False, default=templates.Second())
1071    refDate: datetime.datetime = field(
1072        init=False, repr=False, default=templates.RefDate()
1073    )
1074    productionStatus: templates.Grib2Metadata = field(
1075        init=False, repr=False, default=templates.ProductionStatus()
1076    )
1077    typeOfData: templates.Grib2Metadata = field(
1078        init=False, repr=False, default=templates.TypeOfData()
1079    )
1080
1081    # Section 3 looked up common attributes.  Other looked up attributes are available according
1082    # to the Grid Definition Template.
1083    gridDefinitionSection: NDArray = field(
1084        init=False, repr=False, default=templates.GridDefinitionSection()
1085    )
1086    sourceOfGridDefinition: int = field(
1087        init=False, repr=False, default=templates.SourceOfGridDefinition()
1088    )
1089    numberOfDataPoints: int = field(
1090        init=False, repr=False, default=templates.NumberOfDataPoints()
1091    )
1092    interpretationOfListOfNumbers: templates.Grib2Metadata = field(
1093        init=False, repr=False, default=templates.InterpretationOfListOfNumbers()
1094    )
1095    gridDefinitionTemplateNumber: templates.Grib2Metadata = field(
1096        init=False, repr=False, default=templates.GridDefinitionTemplateNumber()
1097    )
1098    gridDefinitionTemplate: list = field(
1099        init=False, repr=False, default=templates.GridDefinitionTemplate()
1100    )
1101    _earthparams: dict = field(init=False, repr=False, default=templates.EarthParams())
1102    _dxsign: float = field(init=False, repr=False, default=templates.DxSign())
1103    _dysign: float = field(init=False, repr=False, default=templates.DySign())
1104    _llscalefactor: float = field(
1105        init=False, repr=False, default=templates.LLScaleFactor()
1106    )
1107    _lldivisor: float = field(init=False, repr=False, default=templates.LLDivisor())
1108    _xydivisor: float = field(init=False, repr=False, default=templates.XYDivisor())
1109    shapeOfEarth: templates.Grib2Metadata = field(
1110        init=False, repr=False, default=templates.ShapeOfEarth()
1111    )
1112    earthShape: str = field(init=False, repr=False, default=templates.EarthShape())
1113    earthRadius: float = field(init=False, repr=False, default=templates.EarthRadius())
1114    earthMajorAxis: float = field(
1115        init=False, repr=False, default=templates.EarthMajorAxis()
1116    )
1117    earthMinorAxis: float = field(
1118        init=False, repr=False, default=templates.EarthMinorAxis()
1119    )
1120    resolutionAndComponentFlags: list = field(
1121        init=False, repr=False, default=templates.ResolutionAndComponentFlags()
1122    )
1123    ny: int = field(init=False, repr=False, default=templates.Ny())
1124    nx: int = field(init=False, repr=False, default=templates.Nx())
1125    scanModeFlags: list = field(
1126        init=False, repr=False, default=templates.ScanModeFlags()
1127    )
1128    projParameters: dict = field(
1129        init=False, repr=False, default=templates.ProjParameters()
1130    )
1131
1132    # Section 4
1133    productDefinitionTemplateNumber: templates.Grib2Metadata = field(
1134        init=False, repr=False, default=templates.ProductDefinitionTemplateNumber()
1135    )
1136    productDefinitionTemplate: NDArray = field(
1137        init=False, repr=False, default=templates.ProductDefinitionTemplate()
1138    )
1139
1140    # Section 5 looked up common attributes.  Other looked up attributes are
1141    # available according to the Data Representation Template.
1142    numberOfPackedValues: int = field(
1143        init=False, repr=False, default=templates.NumberOfPackedValues()
1144    )
1145    dataRepresentationTemplateNumber: templates.Grib2Metadata = field(
1146        init=False, repr=False, default=templates.DataRepresentationTemplateNumber()
1147    )
1148    dataRepresentationTemplate: list = field(
1149        init=False, repr=False, default=templates.DataRepresentationTemplate()
1150    )
1151    typeOfValues: templates.Grib2Metadata = field(
1152        init=False, repr=False, default=templates.TypeOfValues()
1153    )
1154
1155    def __copy__(self):
1156        """Shallow copy"""
1157        new = Grib2Message(
1158            self.section0,
1159            self.section1,
1160            self.section2,
1161            self.section3,
1162            self.section4,
1163            drtn=self.drtn,
1164        )
1165        return new
1166
1167    def __deepcopy__(self, memo):
1168        """Deep copy"""
1169        new = Grib2Message(
1170            np.copy(self.section0),
1171            np.copy(self.section1),
1172            copy.deepcopy(self.section2),
1173            np.copy(self.section3),
1174            np.copy(self.section4),
1175            np.copy(self.section5),
1176        )
1177        memo[id(self)] = new
1178        new.data = np.copy(self.data)
1179        new.bitmap = None if self.bitmap is None else np.copy(self.bitmap)
1180        return new
1181
1182    def __post_init__(self):
1183        """Set some attributes after init."""
1184        self._auto_nans = _AUTO_NANS
1185        self._coordlist = np.zeros((0), dtype=np.float32)
1186        self._data = None
1187        self._deflist = np.zeros((0), dtype=np.int64)
1188        self._msgnum = -1
1189        self._ondiskarray = None
1190        self._orig_section5 = np.copy(self.section5)
1191        self._signature = self._generate_signature()
1192        try:
1193            self._sha1_section3 = hashlib.sha1(self.section3).hexdigest()
1194        except TypeError:
1195            pass
1196        self.bitMapFlag = templates.Grib2Metadata(self.bitMapFlag, table="6.0")
1197        self.bitmap = None
1198
1199    @property
1200    def _isNDFD(self):
1201        """Check if GRIB2 message is from NWS NDFD"""
1202        return np.all(self.section1[0:2] == [8, 65535])
1203
1204    @property
1205    def _isAerosol(self):
1206        """Check if GRIB2 message contains aerosol data"""
1207        is_aero_template = (
1208            self.productDefinitionTemplateNumber.value in tables.AEROSOL_PDTNS
1209        )
1210        is_aero_param = (self.parameterCategory in {13, 20}) and (
1211            self.parameterNumber in tables.AEROSOL_PARAMS
1212        )
1213        # Check table 4.205 aerosol presence
1214        is_aero_type = self.parameterCategory == 205 and self.parameterNumber == 1
1215        return is_aero_template or is_aero_param or is_aero_type
1216
1217    @property
1218    def _isChemical(self):
1219        """Check if GRIB2 message contains chemical data"""
1220        is_chem_template = (
1221            self.productDefinitionTemplateNumber.value in tables.CHEMICAL_PDTNS
1222        )
1223        is_chem_param = self.parameterCategory == 20
1224        return is_chem_template or is_chem_param
1225
1226    @property
1227    def gdtn(self):
1228        """Return Grid Definition Template Number"""
1229        return self.section3[4]
1230
1231    @property
1232    def gdt(self):
1233        """Return Grid Definition Template."""
1234        return self.gridDefinitionTemplate
1235
1236    @property
1237    def pdtn(self):
1238        """Return Product Definition Template Number."""
1239        return self.section4[1]
1240
1241    @property
1242    def pdt(self):
1243        """Return Product Definition Template."""
1244        return self.productDefinitionTemplate
1245
1246    @property
1247    def drtn(self):
1248        """Return Data Representation Template Number."""
1249        return self.section5[1]
1250
1251    @property
1252    def drt(self):
1253        """Return Data Representation Template."""
1254        return self.dataRepresentationTemplate
1255
1256    @property
1257    def pdy(self):
1258        """Return the PDY ('YYYYMMDD')."""
1259        return "".join([str(i) for i in self.section1[5:8]])
1260
1261    @property
1262    def griddef(self):
1263        """Return a Grib2GridDef instance for a GRIB2 message."""
1264        return Grib2GridDef.from_section3(self.section3)
1265
1266    @property
1267    def lats(self):
1268        """Return grid latitudes."""
1269        return self.latlons()[0]
1270
1271    @property
1272    def lons(self):
1273        """Return grid longitudes."""
1274        return self.latlons()[1]
1275
1276    @property
1277    def min(self):
1278        """Return minimum value of data."""
1279        return np.nanmin(self.data)
1280
1281    @property
1282    def max(self):
1283        """Return maximum value of data."""
1284        return np.nanmax(self.data)
1285
1286    @property
1287    def mean(self):
1288        """Return mean value of data."""
1289        return np.nanmean(self.data)
1290
1291    @property
1292    def median(self):
1293        """Return median value of data."""
1294        return np.nanmedian(self.data)
1295
1296    @property
1297    def shape(self):
1298        """Return shape of data."""
1299        return self.griddef.shape
1300
1301    def __repr__(self):
1302        """
1303        Return an unambiguous string representation of the object.
1304
1305        Returns
1306        -------
1307        repr
1308            A string representation of the object, including information from
1309            sections 0, 1, 3, 4, 5, and 6.
1310        """
1311        info = ""
1312        for sect in [0, 1, 3, 4, 5, 6]:
1313            for k, v in self.attrs_by_section(sect, values=True).items():
1314                info += f"Section {sect}: {k} = {v}\n"
1315        return info
1316
1317    def __str__(self):
1318        """
1319        Return a readable string representation of the object.
1320
1321        Returns
1322        -------
1323        str
1324            A formatted string representation of the object, including
1325            selected attributes.
1326        """
1327        return f"{self._msgnum}:d={self.refDate}:{self.shortName}:{self.fullName} ({self.units}):{self.level}:{self.leadTime}"
1328
1329    def _generate_signature(self):
1330        """Generature SHA-1 hash string from GRIB2 integer sections."""
1331        return hashlib.sha1(
1332            np.concatenate(
1333                (
1334                    self.section0,
1335                    self.section1,
1336                    self.section3,
1337                    self.section4,
1338                    self.section5,
1339                )
1340            )
1341        ).hexdigest()
1342
1343    def attrs_by_section(self, sect: int, values: bool = False):
1344        """
1345        Provide a tuple of attribute names for the given GRIB2 section.
1346
1347        Parameters
1348        ----------
1349        sect
1350            The GRIB2 section number.
1351        values
1352            Optional (default is `False`) argument to return attributes values.
1353
1354        Returns
1355        -------
1356        attrs_by_section
1357            A list of attribute names or dict of name:value pairs if `values =
1358            True`.
1359        """
1360        if sect in {0, 1, 6}:
1361            attrs = templates._section_attrs[sect]
1362        elif sect in {3, 4, 5}:
1363
1364            def _find_class_index(n):
1365                _key = {3: "Grid", 4: "Product", 5: "Data"}
1366                for i, c in enumerate(self.__class__.__mro__):
1367                    if _key[n] in c.__name__:
1368                        return i
1369                else:
1370                    return []
1371
1372            if sys.version_info.minor <= 8:
1373                attrs = templates._section_attrs[sect] + [
1374                    a
1375                    for a in dir(self.__class__.__mro__[_find_class_index(sect)])
1376                    if not a.startswith("_")
1377                ]
1378            else:
1379                attrs = (
1380                    templates._section_attrs[sect]
1381                    + self.__class__.__mro__[_find_class_index(sect)]._attrs()
1382                )
1383        else:
1384            attrs = []
1385        if values:
1386            return {k: getattr(self, k) for k in attrs}
1387        else:
1388            return attrs
1389
1390    def copy(self, deep: bool = True):
1391        """Returns a copy of this Grib2Message.
1392
1393        When `deep=True`, a copy is made of each of the GRIB2 section arrays and
1394        the data are unpacked from the source object and copied into the new
1395        object. Otherwise, a shallow copy of each array is performed and no data
1396        are copied.
1397
1398        Parameters
1399        ----------
1400        deep : bool, default: True
1401            Whether each GRIB2 section array and data are copied onto
1402            the new object. Default is True.
1403
1404        Returns
1405        -------
1406        object : Grib2Message
1407            New Grib2Message object.
1408
1409            .. versionadded:: 2.6.0
1410        """
1411        return copy.deepcopy(self) if deep else copy.copy(self)
1412
1413    def pack(self):
1414        """
1415        Pack GRIB2 section data into a binary message.
1416
1417        It is the user's responsibility to populate the GRIB2 section
1418        information with appropriate metadata.
1419        """
1420        # Create beginning of packed binary message with section 0 and 1 data.
1421        self._sections = []
1422        self._msg, self._pos = g2clib.grib2_create(
1423            self.indicatorSection[2:4], self.identificationSection
1424        )
1425        self._sections += [0, 1]
1426
1427        # Add section 2 if present.
1428        if isinstance(self.section2, bytes) and len(self.section2) > 0:
1429            self._msg, self._pos = g2clib.grib2_addlocal(self._msg, self.section2)
1430            self._sections.append(2)
1431
1432        # Add section 3.
1433        self.section3[1] = self.nx * self.ny
1434        self._msg, self._pos = g2clib.grib2_addgrid(
1435            self._msg,
1436            self.gridDefinitionSection,
1437            self.gridDefinitionTemplate,
1438            self._deflist,
1439        )
1440        self._sections.append(3)
1441
1442        # Prepare data.
1443        if self._data is None:
1444            if self._ondiskarray is None:
1445                raise ValueError(
1446                    "Grib2Message object has no data, thus it cannot be packed."
1447                )
1448        field = np.copy(self.data)
1449        if self.scanModeFlags is not None:
1450            if self.scanModeFlags[3]:
1451                fieldsave = field.astype("f")  # Casting makes a copy
1452                field[1::2, :] = fieldsave[1::2, ::-1]
1453        fld = field.astype("f")
1454
1455        # Prepare bitmap, if necessary
1456        bitmapflag = self.bitMapFlag.value
1457        if bitmapflag == 0:
1458            if self.bitmap is not None:
1459                bmap = np.ravel(self.bitmap).astype(DEFAULT_NUMPY_INT)
1460            else:
1461                bmap = np.ravel(np.where(np.isnan(fld), 0, 1)).astype(DEFAULT_NUMPY_INT)
1462        else:
1463            bmap = None
1464
1465        # Prepare data for packing if nans are present
1466        fld = np.ravel(fld)
1467        if bitmapflag in {0, 254}:
1468            fld = np.where(np.isnan(fld), 0, fld)
1469        else:
1470            if np.isnan(fld).any():
1471                if hasattr(self, "priMissingValue"):
1472                    fld = np.where(np.isnan(fld), self.priMissingValue, fld)
1473            if hasattr(self, "_missvalmap"):
1474                if hasattr(self, "priMissingValue"):
1475                    fld = np.where(self._missvalmap == 1, self.priMissingValue, fld)
1476                if hasattr(self, "secMissingValue"):
1477                    fld = np.where(self._missvalmap == 2, self.secMissingValue, fld)
1478
1479        # Add sections 4, 5, 6, and 7.
1480        self._msg, self._pos = g2clib.grib2_addfield(
1481            self._msg,
1482            self.pdtn,
1483            self.productDefinitionTemplate,
1484            self._coordlist,
1485            self.drtn,
1486            self.dataRepresentationTemplate,
1487            fld,
1488            bitmapflag,
1489            bmap,
1490        )
1491        self._sections.append(4)
1492        self._sections.append(5)
1493        self._sections.append(6)
1494        self._sections.append(7)
1495
1496        # Finalize GRIB2 message with section 8.
1497        self._msg, self._pos = g2clib.grib2_end(self._msg)
1498        self._sections.append(8)
1499        self.section0[-1] = len(self._msg)
1500
1501    @property
1502    def data(self) -> np.array:
1503        """Access the unpacked data values."""
1504        if self._data is None:
1505            if self._auto_nans != _AUTO_NANS:
1506                self._data = self._ondiskarray
1507            self._data = np.asarray(self._ondiskarray)
1508        return self._data
1509
1510    @data.setter
1511    def data(self, arr):
1512        """
1513        Set the internal data array, enforcing shape (ny, nx) and dtype float32.
1514
1515        If the Grid Definition Section (section 3) of Grib2Message object is
1516        not fully formed (i.e. nx, ny = 0, 0), then the shape of the data array
1517        will be used to set nx and ny of the Grib2Message object. It will be the
1518        responsibility of the user to populate the rest of the Grid Definition
1519        Section attributes.
1520
1521        Parameters
1522        ----------
1523        arr : array_like
1524            A 2D array whose shape must match ``(self.ny, self.nx)``.
1525            It will be converted to ``float32`` and C-contiguous if needed.
1526
1527        Raises
1528        ------
1529        ValueError
1530            If the shape of `arr` does not match the expected dimensions.
1531        """
1532        if not isinstance(arr, np.ndarray):
1533            raise ValueError("Grib2Message data only supports numpy arrays")
1534        if self.nx == 0 and self.ny == 0:
1535            self.ny = arr.shape[0]
1536            self.nx = arr.shape[1]
1537        if arr.shape != (self.ny, self.nx):
1538            raise ValueError(
1539                f"Data shape mismatch: expected ({self.ny}, {self.nx}), got {arr.shape}"
1540            )
1541        # Ensure contiguous memory layout (important for C interoperability)
1542        if not arr.flags["C_CONTIGUOUS"]:
1543            arr = np.ascontiguousarray(arr, dtype=np.float32)
1544        self._data = arr
1545
1546    def flush_data(self):
1547        """
1548        Flush the unpacked data values from the Grib2Message object.
1549
1550        Notes
1551        -----
1552        If the Grib2Message object was constructed from "scratch" (i.e.
1553        not read from file), this method will remove the data array from
1554        the object and it cannot be recovered.
1555        """
1556        self._data = None
1557        self.bitmap = None
1558
1559    def __getitem__(self, item):
1560        return self.data[item]
1561
1562    def __setitem__(self, item):
1563        raise NotImplementedError("assignment of data not supported via setitem")
1564
1565    def latlons(self, *args, **kwrgs):
1566        """Alias for `grib2io.Grib2Message.grid` method."""
1567        return self.grid(*args, **kwrgs)
1568
1569    def grid(self, unrotate: bool = True):
1570        """
1571        Return lats,lons (in degrees) of grid.
1572
1573        Currently can handle reg. lat/lon,cglobal Gaussian, mercator,
1574        stereographic, lambert conformal, albers equal-area, space-view and
1575        azimuthal equidistant grids.
1576
1577        Parameters
1578        ----------
1579        unrotate
1580            If `True` [DEFAULT], and grid is rotated lat/lon, then unrotate the
1581            grid, otherwise `False`, do not.
1582
1583        Returns
1584        -------
1585        lats, lons : numpy.ndarray
1586            Returns two numpy.ndarrays with dtype=numpy.float32 of grid
1587            latitudes and longitudes in units of degrees.
1588        """
1589        if self._sha1_section3 in _latlon_datastore.keys():
1590            return (
1591                _latlon_datastore[self._sha1_section3]["latitude"],
1592                _latlon_datastore[self._sha1_section3]["longitude"],
1593            )
1594        gdtn = self.gridDefinitionTemplateNumber.value
1595        reggrid = self.gridDefinitionSection[2] == 0  # This means regular 2-d grid
1596        if gdtn == 0:
1597            # Regular lat/lon grid
1598            lon1, lat1 = self.longitudeFirstGridpoint, self.latitudeFirstGridpoint
1599            lon2, lat2 = self.longitudeLastGridpoint, self.latitudeLastGridpoint
1600            dlon = self.gridlengthXDirection
1601            if lon2 < lon1 and dlon < 0:
1602                lon1 = -lon1
1603            lats = np.linspace(lat1, lat2, self.ny)
1604            if reggrid:
1605                lons = np.linspace(lon1, lon2, self.nx)
1606            else:
1607                lons = np.linspace(lon1, lon2, self.ny * 2)
1608            lons, lats = np.meshgrid(lons, lats)  # Make 2-d arrays.
1609        elif gdtn == 1:  # Rotated Lat/Lon grid
1610            pj = pyproj.Proj(self.projParameters)
1611            lat1, lon1 = self.latitudeFirstGridpoint, self.longitudeFirstGridpoint
1612            lat2, lon2 = self.latitudeLastGridpoint, self.longitudeLastGridpoint
1613            if lon1 > 180.0:
1614                lon1 -= 360.0
1615            if lon2 > 180.0:
1616                lon2 -= 360.0
1617            lats = np.linspace(lat1, lat2, self.ny)
1618            lons = np.linspace(lon1, lon2, self.nx)
1619            lons, lats = np.meshgrid(lons, lats)  # Make 2-d arrays.
1620            if unrotate:
1621                from grib2io.utils import rotated_grid
1622
1623                lats, lons = rotated_grid.unrotate(
1624                    lats,
1625                    lons,
1626                    self.anglePoleRotation,
1627                    self.latitudeSouthernPole,
1628                    self.longitudeSouthernPole,
1629                )
1630        elif gdtn == 40:  # Gaussian grid (only works for global!)
1631            from grib2io.utils.gauss_grid import gaussian_latitudes
1632
1633            lon1, lat1 = self.longitudeFirstGridpoint, self.latitudeFirstGridpoint
1634            lon2, lat2 = self.longitudeLastGridpoint, self.latitudeLastGridpoint
1635            nlats = self.ny
1636            if not reggrid:  # Reduced Gaussian grid.
1637                nlons = 2 * nlats
1638                dlon = 360.0 / nlons
1639            else:
1640                nlons = self.nx
1641                dlon = self.gridlengthXDirection
1642            lons = np.linspace(lon1, lon2, nlons)
1643            # Compute Gaussian lats (north to south)
1644            lats = gaussian_latitudes(nlats)
1645            if lat1 < lat2:  # reverse them if necessary
1646                lats = lats[::-1]
1647            lons, lats = np.meshgrid(lons, lats)
1648        elif gdtn in {10, 20, 30, 31, 110}:
1649            # Mercator, Lambert Conformal, Stereographic, Albers Equal Area,
1650            # Azimuthal Equidistant
1651            dx, dy = self.gridlengthXDirection, self.gridlengthYDirection
1652            lon1, lat1 = self.longitudeFirstGridpoint, self.latitudeFirstGridpoint
1653            pj = pyproj.Proj(self.projParameters)
1654            llcrnrx, llcrnry = pj(lon1, lat1)
1655            x = llcrnrx + dx * np.arange(self.nx)
1656            y = llcrnry + dy * np.arange(self.ny)
1657            x, y = np.meshgrid(x, y)
1658            lons, lats = pj(x, y, inverse=True)
1659        elif gdtn == 90:
1660            # Satellite Projection
1661            dx = self.gridlengthXDirection
1662            dy = self.gridlengthYDirection
1663            pj = pyproj.Proj(self.projParameters)
1664            x = dx * np.indices((self.ny, self.nx), "f")[1, :, :]
1665            x -= 0.5 * x.max()
1666            y = dy * np.indices((self.ny, self.nx), "f")[0, :, :]
1667            y -= 0.5 * y.max()
1668            lons, lats = pj(x, y, inverse=True)
1669            # Set lons,lats to 1.e30 where undefined
1670            abslons = np.fabs(lons)
1671            abslats = np.fabs(lats)
1672            lons = np.where(abslons < 1.0e20, lons, 1.0e30)
1673            lats = np.where(abslats < 1.0e20, lats, 1.0e30)
1674        elif gdtn == 32769:
1675            # Special NCEP Grid, Rotated Lat/Lon, Arakawa E-Grid (Non-Staggered)
1676            from grib2io.utils import arakawa_rotated_grid
1677
1678            di, dj = 0.0, 0.0
1679            do_180 = False
1680            idir = 1 if self.scanModeFlags[0] == 0 else -1
1681            jdir = -1 if self.scanModeFlags[1] == 0 else 1
1682            0 if self.resolutionAndComponentFlags[4] == 0 else 1
1683            la1 = self.latitudeFirstGridpoint
1684            lo1 = self.longitudeFirstGridpoint
1685            clon = self.longitudeCenterGridpoint
1686            clat = self.latitudeCenterGridpoint
1687            lasp = clat - 90.0
1688            losp = clon
1689            llat, llon = arakawa_rotated_grid.ll2rot(la1, lo1, lasp, losp)
1690            la2, lo2 = arakawa_rotated_grid.rot2ll(-llat, -llon, lasp, losp)
1691            rlat = -llat
1692            rlon = -llon
1693            if self.nx == 1:
1694                di = 0.0
1695            elif idir == 1:
1696                ti = rlon
1697                while ti < llon:
1698                    ti += 360.0
1699                di = (ti - llon) / float(self.nx - 1)
1700            else:
1701                ti = llon
1702                while ti < rlon:
1703                    ti += 360.0
1704                di = (ti - rlon) / float(self.nx - 1)
1705            if self.ny == 1:
1706                dj = 0.0
1707            else:
1708                dj = (rlat - llat) / float(self.ny - 1)
1709                if dj < 0.0:
1710                    dj = -dj
1711            if idir == 1:
1712                if llon > rlon:
1713                    llon -= 360.0
1714                if llon < 0 and rlon > 0:
1715                    do_180 = True
1716            else:
1717                if rlon > llon:
1718                    rlon -= 360.0
1719                if rlon < 0 and llon > 0:
1720                    do_180 = True
1721            xlat1d = llat + (np.arange(self.ny) * jdir * dj)
1722            xlon1d = llon + (np.arange(self.nx) * idir * di)
1723            xlons, xlats = np.meshgrid(xlon1d, xlat1d)
1724            rot2ll_vectorized = np.vectorize(arakawa_rotated_grid.rot2ll)
1725            lats, lons = rot2ll_vectorized(xlats, xlons, lasp, losp)
1726            if do_180:
1727                lons = np.where(lons > 180.0, lons - 360.0, lons)
1728            vector_rotation_angles_vectorized = np.vectorize(
1729                arakawa_rotated_grid.vector_rotation_angles
1730            )
1731            rots = vector_rotation_angles_vectorized(lats, lons, clat, losp, xlats)
1732            del xlat1d, xlon1d, xlats, xlons
1733        else:
1734            raise ValueError("Unsupported grid")
1735
1736        _latlon_datastore[self._sha1_section3] = dict(latitude=lats, longitude=lons)
1737        try:
1738            _latlon_datastore[self._sha1_section3]["vector_rotation_angles"] = rots
1739        except NameError:
1740            pass
1741
1742        return lats, lons
1743
1744    def map_keys(self):
1745        """
1746        Unpack data grid replacing integer values with strings.
1747
1748        These types of fields are categorical or classifications where data
1749        values do not represent an observable or predictable physical quantity.
1750        An example of such a field would be [Dominant Precipitation Type -
1751        DPTYPE](https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_table4-201.shtml)
1752
1753        Returns
1754        -------
1755        map_keys
1756            numpy.ndarray of string values per element.
1757        """
1758        hold_auto_nans = _AUTO_NANS
1759        set_auto_nans(False)
1760
1761        if (np.all(self.section1[0:2] == [7, 14]) and self.shortName == "PWTHER") or (
1762            self._isNDFD and self.shortName in {"WX", "WWA"}
1763        ):
1764            keys = utils.decode_wx_strings(self.section2)
1765            if hasattr(self, "priMissingValue") and self.priMissingValue not in [
1766                None,
1767                0,
1768            ]:
1769                keys[int(self.priMissingValue)] = "Missing"
1770            if hasattr(self, "secMissingValue") and self.secMissingValue not in [
1771                None,
1772                0,
1773            ]:
1774                keys[int(self.secMissingValue)] = "Missing"
1775            u, inv = np.unique(self.data, return_inverse=True)
1776            fld = np.array([keys[x] for x in u])[inv].reshape(self.data.shape)
1777        else:
1778            # For data whose units are defined in a code table (i.e. classification or mask)
1779            tblname = re.findall(r"\d\.\d+", self.units, re.IGNORECASE)[0]
1780            fld = self.data.astype(np.int32).astype(str)
1781            tbl = tables.get_table(tblname, expand=True)
1782            for val in np.unique(fld):
1783                fld = np.where(fld == val, tbl[val], fld)
1784        set_auto_nans(hold_auto_nans)
1785        return fld
1786
1787    def to_bytes(self, validate: bool = True):
1788        """
1789        Return packed GRIB2 message in bytes format.
1790
1791        This will be useful for exporting data in non-file formats. For example,
1792        can be used to output grib data directly to S3 using the boto3 client
1793        without the need to write a temporary file to upload first.
1794
1795        Parameters
1796        ----------
1797        validate: default=True
1798            If `True` (DEFAULT), validates first/last four bytes for proper
1799            formatting, else returns None. If `False`, message is output as is.
1800
1801        Returns
1802        -------
1803        to_bytes
1804            Returns GRIB2 formatted message as bytes.
1805        """
1806        if hasattr(self, "_msg"):
1807            if validate:
1808                if self.validate():
1809                    return self._msg
1810                else:
1811                    return None
1812            else:
1813                return self._msg
1814        else:
1815            return None
1816
1817    def interpolate(
1818        self, method, grid_def_out, method_options=None, drtn=None, num_threads=1
1819    ):
1820        """
1821        Grib2Message Interpolator
1822
1823        Performs spatial interpolation via the [NCEPLIBS-ip
1824        library](https://github.com/NOAA-EMC/NCEPLIBS-ip). This interpolate
1825        method only supports scalar interpolation. If you need to perform
1826        vector interpolation, use the module-level `grib2io.interpolate`
1827        function.
1828
1829        Parameters
1830        ----------
1831        method
1832            Interpolate method to use. This can either be an integer or string
1833            using the following mapping:
1834
1835            | Interpolate Scheme | Integer Value |
1836            | :---:              | :---:         |
1837            | 'bilinear'         | 0             |
1838            | 'bicubic'          | 1             |
1839            | 'neighbor'         | 2             |
1840            | 'budget'           | 3             |
1841            | 'spectral'         | 4             |
1842            | 'neighbor-budget'  | 6             |
1843
1844        grid_def_out : grib2io.Grib2GridDef
1845            Grib2GridDef object of the output grid.
1846        method_options : list of ints, optional
1847            Interpolation options. See the NCEPLIBS-ip documentation for
1848            more information on how these are used.
1849        drtn
1850            Data Representation Template to be used for the returned
1851            interpolated GRIB2 message. When `None`, the data representation
1852            template of the source GRIB2 message is used. Once again, it is the
1853            user's responsibility to properly set the Data Representation
1854            Template attributes.
1855        num_threads : int, optional
1856            Number of OpenMP threads to use for interpolation. The default
1857            value is 1. If NCEPLIBS-ip and grib2io's iplib extension module
1858            was not built with OpenMP, then this keyword argument and value
1859            will have no impact.
1860
1861        Returns
1862        -------
1863        interpolate
1864            If interpolating to a grid, a new Grib2Message object is returned.
1865            The GRIB2 metadata of the new Grib2Message object is identical to
1866            the input except where required to be different because of the new
1867            grid specs and possibly a new data representation template.
1868
1869            If interpolating to station points, the interpolated data values are
1870            returned as a numpy.ndarray.
1871        """
1872        section0 = self.section0
1873        section0[-1] = 0
1874        gds = [0, grid_def_out.npoints, 0, 255, grid_def_out.gdtn]
1875        section3 = np.concatenate((gds, grid_def_out.gdt))
1876        drtn = self.drtn if drtn is None else drtn
1877
1878        msg = Grib2Message(
1879            section0,
1880            self.section1,
1881            self.section2,
1882            section3,
1883            self.section4,
1884            None,
1885            self.bitMapFlag.value,
1886            drtn=drtn,
1887        )
1888
1889        msg._msgnum = -1
1890        msg._deflist = self._deflist
1891        msg._coordlist = self._coordlist
1892        if msg.typeOfValues == 0:
1893            pass
1894        elif msg.typeOfValues == 1:
1895            pass
1896        msg._data = interpolate(
1897            self.data,
1898            method,
1899            Grib2GridDef.from_section3(self.section3),
1900            grid_def_out,
1901            method_options=method_options,
1902            num_threads=num_threads,
1903        ).reshape(msg.ny, msg.nx)
1904        msg.section5[0] = grid_def_out.npoints
1905        return msg
1906
1907    def subset(self, lats, lons):
1908        """
1909        Return a spatial subset.
1910
1911        Currently only supports regular grids of the following types:
1912
1913        | Grid Type                                                    | gdtn  |
1914        | :---:                                                        | :---: |
1915        | Latitude/Longitude, Equidistant Cylindrical, or Plate Carree | 0     |
1916        | Rotated Latitude/Longitude                                   | 1     |
1917        | Mercator                                                     | 10    |
1918        | Polar Stereographic                                          | 20    |
1919        | Lambert Conformal                                            | 30    |
1920        | Albers Equal-Area                                            | 31    |
1921        | Gaussian Latitude/Longitude                                  | 40    |
1922        | Equatorial Azimuthal Equidistant Projection                  | 110   |
1923
1924        Parameters
1925        ----------
1926        lats
1927            List or tuple of latitudes.  The minimum and maximum latitudes will
1928            be used to define the southern and northern boundaries.
1929
1930            The order of the latitudes is not important.  The function will
1931            determine which is the minimum and maximum.
1932
1933            The latitudes should be in decimal degrees with 0.0 at the equator,
1934            positive values in the northern hemisphere increasing to 90, and
1935            negative values in the southern hemisphere decreasing to -90.
1936        lons
1937            List or tuple of longitudes.  The minimum and maximum longitudes
1938            will be used to define the western and eastern boundaries.
1939
1940            The order of the longitudes is not important.  The function will
1941            determine which is the minimum and maximum.
1942
1943            GRIB2 longitudes should be in decimal degrees with 0.0 at the prime
1944            meridian, positive values increasing eastward to 360.  There are no
1945            negative GRIB2 longitudes.
1946
1947            The typical west longitudes that start at 0.0 at the prime meridian
1948            and decrease to -180 westward, are converted to GRIB2 longitudes by
1949            '360 - (absolute value of the west longitude)' where typical
1950            eastern longitudes are unchanged as GRIB2 longitudes.
1951
1952        Returns
1953        -------
1954        subset
1955            A spatial subset of a GRIB2 message.
1956        """
1957        if self.gdtn not in [0, 1, 10, 20, 30, 31, 40, 110]:
1958            raise ValueError(
1959                """
1960
1961Subset only works for
1962    Latitude/Longitude, Equidistant Cylindrical, or Plate Carree (gdtn=0)
1963    Rotated Latitude/Longitude (gdtn=1)
1964    Mercator (gdtn=10)
1965    Polar Stereographic (gdtn=20)
1966    Lambert Conformal (gdtn=30)
1967    Albers Equal-Area (gdtn=31)
1968    Gaussian Latitude/Longitude (gdtn=40)
1969    Equatorial Azimuthal Equidistant Projection (gdtn=110)
1970
1971"""
1972            )
1973
1974        if self.nx == 0 or self.ny == 0:
1975            raise ValueError(
1976                """
1977
1978Subset only works for regular grids.
1979
1980"""
1981            )
1982
1983        newmsg = Grib2Message(
1984            np.copy(self.section0),
1985            np.copy(self.section1),
1986            np.copy(self.section2),
1987            np.copy(self.section3),
1988            np.copy(self.section4),
1989            np.copy(self.section5),
1990        )
1991
1992        msglats, msglons = self.grid()
1993
1994        la1 = np.max(lats)
1995        lo1 = np.min(lons)
1996        la2 = np.min(lats)
1997        lo2 = np.max(lons)
1998
1999        # Find the indices of the first and last grid points to the nearest
2000        # lat/lon values in the grid.
2001        first_lat = np.abs(msglats - la1)
2002        first_lon = np.abs(msglons - lo1)
2003        max_idx = np.maximum(first_lat, first_lon)
2004        first_j, first_i = np.where(max_idx == np.min(max_idx))
2005
2006        last_lat = np.abs(msglats - la2)
2007        last_lon = np.abs(msglons - lo2)
2008        max_idx = np.maximum(last_lat, last_lon)
2009        last_j, last_i = np.where(max_idx == np.min(max_idx))
2010
2011        setattr(newmsg, "latitudeFirstGridpoint", msglats[first_j[0], first_i[0]])
2012        setattr(newmsg, "longitudeFirstGridpoint", msglons[first_j[0], first_i[0]])
2013        setattr(newmsg, "nx", np.abs(first_i[0] - last_i[0]))
2014        setattr(newmsg, "ny", np.abs(first_j[0] - last_j[0]))
2015
2016        # Set *LastGridpoint attributes even if only used for gdtn=[0, 1, 40].
2017        # This information is used to subset xarray datasets and even though
2018        # unnecessary for some supported grid types, it won't affect a grib2io
2019        # message to set them.
2020        setattr(newmsg, "latitudeLastGridpoint", msglats[last_j[0], last_i[0]])
2021        setattr(newmsg, "longitudeLastGridpoint", msglons[last_j[0], last_i[0]])
2022
2023        setattr(
2024            newmsg,
2025            "data",
2026            self.data[
2027                min(first_j[0], last_j[0]) : max(first_j[0], last_j[0]),
2028                min(first_i[0], last_i[0]) : max(first_i[0], last_i[0]),
2029            ].copy(),
2030        )
2031
2032        # Need to set the newmsg._sha1_section3 to a blank string so the grid
2033        # method ignores the cached lat/lon values.  This will force the grid
2034        # method to recompute the lat/lon values for the subsetted grid.
2035        newmsg._sha1_section3 = ""
2036        newmsg.grid()
2037
2038        return newmsg
2039
2040    def validate(self):
2041        """
2042        Validate a complete GRIB2 message.
2043
2044        The g2c library does its own internal validation when g2_gribend() is called, but
2045        we will check in grib2io also. The validation checks if the first 4 bytes in
2046        self._msg is 'GRIB' and '7777' as the last 4 bytes and that the message length in
2047        section 0 equals the length of the packed message.
2048
2049        Returns
2050        -------
2051        `True` if the packed GRIB2 message is complete and well-formed, `False` otherwise.
2052        """
2053        valid = False
2054        if hasattr(self, "_msg"):
2055            if self._msg[0:4] + self._msg[-4:] == b"GRIB7777":
2056                if self.section0[-1] == len(self._msg):
2057                    valid = True
2058        return valid

GRIB2 Message base class.

_Grib2Message( section0: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]], section1: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]], section2: bytes, section3: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]], section4: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]], section5: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]], bitMapFlag: grib2io.templates.Grib2Metadata = 255)
section0: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]]
section1: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]]
section2: bytes
section3: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]]
section4: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]]
section5: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]]
indicatorSection: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]]
identificationSection: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]]

GRIB2 Section 1, Identification Section

year: int

Year of reference time

month: int

Month of reference time

day: int

Day of reference time

hour: int

Hour of reference time

minute: int

Minute of reference time

second: int

Second of reference time

refDate: datetime.datetime

Reference Date. NOTE: This is a datetime.datetime object.

gridDefinitionSection: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]]

GRIB2 Section 3, Grid Definition Section

sourceOfGridDefinition: int
numberOfDataPoints: int

Number of Data Points

interpretationOfListOfNumbers: grib2io.templates.Grib2Metadata

Interpretation of List of Numbers

gridDefinitionTemplate: list

Grid definition template

earthShape: str

Description of the shape of the Earth

earthRadius: float

Radius of the Earth (Assumes "spherical")

earthMajorAxis: float

Major Axis of the Earth (Assumes "oblate spheroid" or "ellipsoid")

earthMinorAxis: float

Minor Axis of the Earth (Assumes "oblate spheroid" or "ellipsoid")

resolutionAndComponentFlags: list
ny: int

Number of grid points in the Y-direction (generally North-South)

nx: int

Number of grid points in the X-direction (generally East-West)

scanModeFlags: list
projParameters: dict

PROJ Parameters to define the reference system

productDefinitionTemplate: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]]

Product Definition Template

numberOfPackedValues: int

Number of Packed Values

dataRepresentationTemplate: list

Data Representation Template

gdtn
1226    @property
1227    def gdtn(self):
1228        """Return Grid Definition Template Number"""
1229        return self.section3[4]

Return Grid Definition Template Number

gdt
1231    @property
1232    def gdt(self):
1233        """Return Grid Definition Template."""
1234        return self.gridDefinitionTemplate

Return Grid Definition Template.

pdtn
1236    @property
1237    def pdtn(self):
1238        """Return Product Definition Template Number."""
1239        return self.section4[1]

Return Product Definition Template Number.

pdt
1241    @property
1242    def pdt(self):
1243        """Return Product Definition Template."""
1244        return self.productDefinitionTemplate

Return Product Definition Template.

drtn
1246    @property
1247    def drtn(self):
1248        """Return Data Representation Template Number."""
1249        return self.section5[1]

Return Data Representation Template Number.

drt
1251    @property
1252    def drt(self):
1253        """Return Data Representation Template."""
1254        return self.dataRepresentationTemplate

Return Data Representation Template.

pdy
1256    @property
1257    def pdy(self):
1258        """Return the PDY ('YYYYMMDD')."""
1259        return "".join([str(i) for i in self.section1[5:8]])

Return the PDY ('YYYYMMDD').

griddef
1261    @property
1262    def griddef(self):
1263        """Return a Grib2GridDef instance for a GRIB2 message."""
1264        return Grib2GridDef.from_section3(self.section3)

Return a Grib2GridDef instance for a GRIB2 message.

lats
1266    @property
1267    def lats(self):
1268        """Return grid latitudes."""
1269        return self.latlons()[0]

Return grid latitudes.

lons
1271    @property
1272    def lons(self):
1273        """Return grid longitudes."""
1274        return self.latlons()[1]

Return grid longitudes.

min
1276    @property
1277    def min(self):
1278        """Return minimum value of data."""
1279        return np.nanmin(self.data)

Return minimum value of data.

max
1281    @property
1282    def max(self):
1283        """Return maximum value of data."""
1284        return np.nanmax(self.data)

Return maximum value of data.

mean
1286    @property
1287    def mean(self):
1288        """Return mean value of data."""
1289        return np.nanmean(self.data)

Return mean value of data.

median
1291    @property
1292    def median(self):
1293        """Return median value of data."""
1294        return np.nanmedian(self.data)

Return median value of data.

shape
1296    @property
1297    def shape(self):
1298        """Return shape of data."""
1299        return self.griddef.shape

Return shape of data.

def attrs_by_section(self, sect: int, values: bool = False):
1343    def attrs_by_section(self, sect: int, values: bool = False):
1344        """
1345        Provide a tuple of attribute names for the given GRIB2 section.
1346
1347        Parameters
1348        ----------
1349        sect
1350            The GRIB2 section number.
1351        values
1352            Optional (default is `False`) argument to return attributes values.
1353
1354        Returns
1355        -------
1356        attrs_by_section
1357            A list of attribute names or dict of name:value pairs if `values =
1358            True`.
1359        """
1360        if sect in {0, 1, 6}:
1361            attrs = templates._section_attrs[sect]
1362        elif sect in {3, 4, 5}:
1363
1364            def _find_class_index(n):
1365                _key = {3: "Grid", 4: "Product", 5: "Data"}
1366                for i, c in enumerate(self.__class__.__mro__):
1367                    if _key[n] in c.__name__:
1368                        return i
1369                else:
1370                    return []
1371
1372            if sys.version_info.minor <= 8:
1373                attrs = templates._section_attrs[sect] + [
1374                    a
1375                    for a in dir(self.__class__.__mro__[_find_class_index(sect)])
1376                    if not a.startswith("_")
1377                ]
1378            else:
1379                attrs = (
1380                    templates._section_attrs[sect]
1381                    + self.__class__.__mro__[_find_class_index(sect)]._attrs()
1382                )
1383        else:
1384            attrs = []
1385        if values:
1386            return {k: getattr(self, k) for k in attrs}
1387        else:
1388            return attrs

Provide a tuple of attribute names for the given GRIB2 section.

Parameters
  • sect: The GRIB2 section number.
  • values: Optional (default is False) argument to return attributes values.
Returns
  • attrs_by_section: A list of attribute names or dict of name:value pairs if values = True.
def copy(self, deep: bool = True):
1390    def copy(self, deep: bool = True):
1391        """Returns a copy of this Grib2Message.
1392
1393        When `deep=True`, a copy is made of each of the GRIB2 section arrays and
1394        the data are unpacked from the source object and copied into the new
1395        object. Otherwise, a shallow copy of each array is performed and no data
1396        are copied.
1397
1398        Parameters
1399        ----------
1400        deep : bool, default: True
1401            Whether each GRIB2 section array and data are copied onto
1402            the new object. Default is True.
1403
1404        Returns
1405        -------
1406        object : Grib2Message
1407            New Grib2Message object.
1408
1409            .. versionadded:: 2.6.0
1410        """
1411        return copy.deepcopy(self) if deep else copy.copy(self)

Returns a copy of this Grib2Message.

When deep=True, a copy is made of each of the GRIB2 section arrays and the data are unpacked from the source object and copied into the new object. Otherwise, a shallow copy of each array is performed and no data are copied.

Parameters
  • deep : bool, default (True): Whether each GRIB2 section array and data are copied onto the new object. Default is True.
Returns
  • object (Grib2Message): New Grib2Message object.

    New in version 2.6.0.

def pack(self):
1413    def pack(self):
1414        """
1415        Pack GRIB2 section data into a binary message.
1416
1417        It is the user's responsibility to populate the GRIB2 section
1418        information with appropriate metadata.
1419        """
1420        # Create beginning of packed binary message with section 0 and 1 data.
1421        self._sections = []
1422        self._msg, self._pos = g2clib.grib2_create(
1423            self.indicatorSection[2:4], self.identificationSection
1424        )
1425        self._sections += [0, 1]
1426
1427        # Add section 2 if present.
1428        if isinstance(self.section2, bytes) and len(self.section2) > 0:
1429            self._msg, self._pos = g2clib.grib2_addlocal(self._msg, self.section2)
1430            self._sections.append(2)
1431
1432        # Add section 3.
1433        self.section3[1] = self.nx * self.ny
1434        self._msg, self._pos = g2clib.grib2_addgrid(
1435            self._msg,
1436            self.gridDefinitionSection,
1437            self.gridDefinitionTemplate,
1438            self._deflist,
1439        )
1440        self._sections.append(3)
1441
1442        # Prepare data.
1443        if self._data is None:
1444            if self._ondiskarray is None:
1445                raise ValueError(
1446                    "Grib2Message object has no data, thus it cannot be packed."
1447                )
1448        field = np.copy(self.data)
1449        if self.scanModeFlags is not None:
1450            if self.scanModeFlags[3]:
1451                fieldsave = field.astype("f")  # Casting makes a copy
1452                field[1::2, :] = fieldsave[1::2, ::-1]
1453        fld = field.astype("f")
1454
1455        # Prepare bitmap, if necessary
1456        bitmapflag = self.bitMapFlag.value
1457        if bitmapflag == 0:
1458            if self.bitmap is not None:
1459                bmap = np.ravel(self.bitmap).astype(DEFAULT_NUMPY_INT)
1460            else:
1461                bmap = np.ravel(np.where(np.isnan(fld), 0, 1)).astype(DEFAULT_NUMPY_INT)
1462        else:
1463            bmap = None
1464
1465        # Prepare data for packing if nans are present
1466        fld = np.ravel(fld)
1467        if bitmapflag in {0, 254}:
1468            fld = np.where(np.isnan(fld), 0, fld)
1469        else:
1470            if np.isnan(fld).any():
1471                if hasattr(self, "priMissingValue"):
1472                    fld = np.where(np.isnan(fld), self.priMissingValue, fld)
1473            if hasattr(self, "_missvalmap"):
1474                if hasattr(self, "priMissingValue"):
1475                    fld = np.where(self._missvalmap == 1, self.priMissingValue, fld)
1476                if hasattr(self, "secMissingValue"):
1477                    fld = np.where(self._missvalmap == 2, self.secMissingValue, fld)
1478
1479        # Add sections 4, 5, 6, and 7.
1480        self._msg, self._pos = g2clib.grib2_addfield(
1481            self._msg,
1482            self.pdtn,
1483            self.productDefinitionTemplate,
1484            self._coordlist,
1485            self.drtn,
1486            self.dataRepresentationTemplate,
1487            fld,
1488            bitmapflag,
1489            bmap,
1490        )
1491        self._sections.append(4)
1492        self._sections.append(5)
1493        self._sections.append(6)
1494        self._sections.append(7)
1495
1496        # Finalize GRIB2 message with section 8.
1497        self._msg, self._pos = g2clib.grib2_end(self._msg)
1498        self._sections.append(8)
1499        self.section0[-1] = len(self._msg)

Pack GRIB2 section data into a binary message.

It is the user's responsibility to populate the GRIB2 section information with appropriate metadata.

data: <built-in function array>
1501    @property
1502    def data(self) -> np.array:
1503        """Access the unpacked data values."""
1504        if self._data is None:
1505            if self._auto_nans != _AUTO_NANS:
1506                self._data = self._ondiskarray
1507            self._data = np.asarray(self._ondiskarray)
1508        return self._data

Access the unpacked data values.

def flush_data(self):
1546    def flush_data(self):
1547        """
1548        Flush the unpacked data values from the Grib2Message object.
1549
1550        Notes
1551        -----
1552        If the Grib2Message object was constructed from "scratch" (i.e.
1553        not read from file), this method will remove the data array from
1554        the object and it cannot be recovered.
1555        """
1556        self._data = None
1557        self.bitmap = None

Flush the unpacked data values from the Grib2Message object.

Notes

If the Grib2Message object was constructed from "scratch" (i.e. not read from file), this method will remove the data array from the object and it cannot be recovered.

def latlons(self, *args, **kwrgs):
1565    def latlons(self, *args, **kwrgs):
1566        """Alias for `grib2io.Grib2Message.grid` method."""
1567        return self.grid(*args, **kwrgs)

Alias for grib2io.Grib2Message.grid method.

def grid(self, unrotate: bool = True):
1569    def grid(self, unrotate: bool = True):
1570        """
1571        Return lats,lons (in degrees) of grid.
1572
1573        Currently can handle reg. lat/lon,cglobal Gaussian, mercator,
1574        stereographic, lambert conformal, albers equal-area, space-view and
1575        azimuthal equidistant grids.
1576
1577        Parameters
1578        ----------
1579        unrotate
1580            If `True` [DEFAULT], and grid is rotated lat/lon, then unrotate the
1581            grid, otherwise `False`, do not.
1582
1583        Returns
1584        -------
1585        lats, lons : numpy.ndarray
1586            Returns two numpy.ndarrays with dtype=numpy.float32 of grid
1587            latitudes and longitudes in units of degrees.
1588        """
1589        if self._sha1_section3 in _latlon_datastore.keys():
1590            return (
1591                _latlon_datastore[self._sha1_section3]["latitude"],
1592                _latlon_datastore[self._sha1_section3]["longitude"],
1593            )
1594        gdtn = self.gridDefinitionTemplateNumber.value
1595        reggrid = self.gridDefinitionSection[2] == 0  # This means regular 2-d grid
1596        if gdtn == 0:
1597            # Regular lat/lon grid
1598            lon1, lat1 = self.longitudeFirstGridpoint, self.latitudeFirstGridpoint
1599            lon2, lat2 = self.longitudeLastGridpoint, self.latitudeLastGridpoint
1600            dlon = self.gridlengthXDirection
1601            if lon2 < lon1 and dlon < 0:
1602                lon1 = -lon1
1603            lats = np.linspace(lat1, lat2, self.ny)
1604            if reggrid:
1605                lons = np.linspace(lon1, lon2, self.nx)
1606            else:
1607                lons = np.linspace(lon1, lon2, self.ny * 2)
1608            lons, lats = np.meshgrid(lons, lats)  # Make 2-d arrays.
1609        elif gdtn == 1:  # Rotated Lat/Lon grid
1610            pj = pyproj.Proj(self.projParameters)
1611            lat1, lon1 = self.latitudeFirstGridpoint, self.longitudeFirstGridpoint
1612            lat2, lon2 = self.latitudeLastGridpoint, self.longitudeLastGridpoint
1613            if lon1 > 180.0:
1614                lon1 -= 360.0
1615            if lon2 > 180.0:
1616                lon2 -= 360.0
1617            lats = np.linspace(lat1, lat2, self.ny)
1618            lons = np.linspace(lon1, lon2, self.nx)
1619            lons, lats = np.meshgrid(lons, lats)  # Make 2-d arrays.
1620            if unrotate:
1621                from grib2io.utils import rotated_grid
1622
1623                lats, lons = rotated_grid.unrotate(
1624                    lats,
1625                    lons,
1626                    self.anglePoleRotation,
1627                    self.latitudeSouthernPole,
1628                    self.longitudeSouthernPole,
1629                )
1630        elif gdtn == 40:  # Gaussian grid (only works for global!)
1631            from grib2io.utils.gauss_grid import gaussian_latitudes
1632
1633            lon1, lat1 = self.longitudeFirstGridpoint, self.latitudeFirstGridpoint
1634            lon2, lat2 = self.longitudeLastGridpoint, self.latitudeLastGridpoint
1635            nlats = self.ny
1636            if not reggrid:  # Reduced Gaussian grid.
1637                nlons = 2 * nlats
1638                dlon = 360.0 / nlons
1639            else:
1640                nlons = self.nx
1641                dlon = self.gridlengthXDirection
1642            lons = np.linspace(lon1, lon2, nlons)
1643            # Compute Gaussian lats (north to south)
1644            lats = gaussian_latitudes(nlats)
1645            if lat1 < lat2:  # reverse them if necessary
1646                lats = lats[::-1]
1647            lons, lats = np.meshgrid(lons, lats)
1648        elif gdtn in {10, 20, 30, 31, 110}:
1649            # Mercator, Lambert Conformal, Stereographic, Albers Equal Area,
1650            # Azimuthal Equidistant
1651            dx, dy = self.gridlengthXDirection, self.gridlengthYDirection
1652            lon1, lat1 = self.longitudeFirstGridpoint, self.latitudeFirstGridpoint
1653            pj = pyproj.Proj(self.projParameters)
1654            llcrnrx, llcrnry = pj(lon1, lat1)
1655            x = llcrnrx + dx * np.arange(self.nx)
1656            y = llcrnry + dy * np.arange(self.ny)
1657            x, y = np.meshgrid(x, y)
1658            lons, lats = pj(x, y, inverse=True)
1659        elif gdtn == 90:
1660            # Satellite Projection
1661            dx = self.gridlengthXDirection
1662            dy = self.gridlengthYDirection
1663            pj = pyproj.Proj(self.projParameters)
1664            x = dx * np.indices((self.ny, self.nx), "f")[1, :, :]
1665            x -= 0.5 * x.max()
1666            y = dy * np.indices((self.ny, self.nx), "f")[0, :, :]
1667            y -= 0.5 * y.max()
1668            lons, lats = pj(x, y, inverse=True)
1669            # Set lons,lats to 1.e30 where undefined
1670            abslons = np.fabs(lons)
1671            abslats = np.fabs(lats)
1672            lons = np.where(abslons < 1.0e20, lons, 1.0e30)
1673            lats = np.where(abslats < 1.0e20, lats, 1.0e30)
1674        elif gdtn == 32769:
1675            # Special NCEP Grid, Rotated Lat/Lon, Arakawa E-Grid (Non-Staggered)
1676            from grib2io.utils import arakawa_rotated_grid
1677
1678            di, dj = 0.0, 0.0
1679            do_180 = False
1680            idir = 1 if self.scanModeFlags[0] == 0 else -1
1681            jdir = -1 if self.scanModeFlags[1] == 0 else 1
1682            0 if self.resolutionAndComponentFlags[4] == 0 else 1
1683            la1 = self.latitudeFirstGridpoint
1684            lo1 = self.longitudeFirstGridpoint
1685            clon = self.longitudeCenterGridpoint
1686            clat = self.latitudeCenterGridpoint
1687            lasp = clat - 90.0
1688            losp = clon
1689            llat, llon = arakawa_rotated_grid.ll2rot(la1, lo1, lasp, losp)
1690            la2, lo2 = arakawa_rotated_grid.rot2ll(-llat, -llon, lasp, losp)
1691            rlat = -llat
1692            rlon = -llon
1693            if self.nx == 1:
1694                di = 0.0
1695            elif idir == 1:
1696                ti = rlon
1697                while ti < llon:
1698                    ti += 360.0
1699                di = (ti - llon) / float(self.nx - 1)
1700            else:
1701                ti = llon
1702                while ti < rlon:
1703                    ti += 360.0
1704                di = (ti - rlon) / float(self.nx - 1)
1705            if self.ny == 1:
1706                dj = 0.0
1707            else:
1708                dj = (rlat - llat) / float(self.ny - 1)
1709                if dj < 0.0:
1710                    dj = -dj
1711            if idir == 1:
1712                if llon > rlon:
1713                    llon -= 360.0
1714                if llon < 0 and rlon > 0:
1715                    do_180 = True
1716            else:
1717                if rlon > llon:
1718                    rlon -= 360.0
1719                if rlon < 0 and llon > 0:
1720                    do_180 = True
1721            xlat1d = llat + (np.arange(self.ny) * jdir * dj)
1722            xlon1d = llon + (np.arange(self.nx) * idir * di)
1723            xlons, xlats = np.meshgrid(xlon1d, xlat1d)
1724            rot2ll_vectorized = np.vectorize(arakawa_rotated_grid.rot2ll)
1725            lats, lons = rot2ll_vectorized(xlats, xlons, lasp, losp)
1726            if do_180:
1727                lons = np.where(lons > 180.0, lons - 360.0, lons)
1728            vector_rotation_angles_vectorized = np.vectorize(
1729                arakawa_rotated_grid.vector_rotation_angles
1730            )
1731            rots = vector_rotation_angles_vectorized(lats, lons, clat, losp, xlats)
1732            del xlat1d, xlon1d, xlats, xlons
1733        else:
1734            raise ValueError("Unsupported grid")
1735
1736        _latlon_datastore[self._sha1_section3] = dict(latitude=lats, longitude=lons)
1737        try:
1738            _latlon_datastore[self._sha1_section3]["vector_rotation_angles"] = rots
1739        except NameError:
1740            pass
1741
1742        return lats, lons

Return lats,lons (in degrees) of grid.

Currently can handle reg. lat/lon,cglobal Gaussian, mercator, stereographic, lambert conformal, albers equal-area, space-view and azimuthal equidistant grids.

Parameters
  • unrotate: If True [DEFAULT], and grid is rotated lat/lon, then unrotate the grid, otherwise False, do not.
Returns
  • lats, lons (numpy.ndarray): Returns two numpy.ndarrays with dtype=numpy.float32 of grid latitudes and longitudes in units of degrees.
def map_keys(self):
1744    def map_keys(self):
1745        """
1746        Unpack data grid replacing integer values with strings.
1747
1748        These types of fields are categorical or classifications where data
1749        values do not represent an observable or predictable physical quantity.
1750        An example of such a field would be [Dominant Precipitation Type -
1751        DPTYPE](https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_table4-201.shtml)
1752
1753        Returns
1754        -------
1755        map_keys
1756            numpy.ndarray of string values per element.
1757        """
1758        hold_auto_nans = _AUTO_NANS
1759        set_auto_nans(False)
1760
1761        if (np.all(self.section1[0:2] == [7, 14]) and self.shortName == "PWTHER") or (
1762            self._isNDFD and self.shortName in {"WX", "WWA"}
1763        ):
1764            keys = utils.decode_wx_strings(self.section2)
1765            if hasattr(self, "priMissingValue") and self.priMissingValue not in [
1766                None,
1767                0,
1768            ]:
1769                keys[int(self.priMissingValue)] = "Missing"
1770            if hasattr(self, "secMissingValue") and self.secMissingValue not in [
1771                None,
1772                0,
1773            ]:
1774                keys[int(self.secMissingValue)] = "Missing"
1775            u, inv = np.unique(self.data, return_inverse=True)
1776            fld = np.array([keys[x] for x in u])[inv].reshape(self.data.shape)
1777        else:
1778            # For data whose units are defined in a code table (i.e. classification or mask)
1779            tblname = re.findall(r"\d\.\d+", self.units, re.IGNORECASE)[0]
1780            fld = self.data.astype(np.int32).astype(str)
1781            tbl = tables.get_table(tblname, expand=True)
1782            for val in np.unique(fld):
1783                fld = np.where(fld == val, tbl[val], fld)
1784        set_auto_nans(hold_auto_nans)
1785        return fld

Unpack data grid replacing integer values with strings.

These types of fields are categorical or classifications where data values do not represent an observable or predictable physical quantity. An example of such a field would be Dominant Precipitation Type - DPTYPE

Returns
  • map_keys: numpy.ndarray of string values per element.
def to_bytes(self, validate: bool = True):
1787    def to_bytes(self, validate: bool = True):
1788        """
1789        Return packed GRIB2 message in bytes format.
1790
1791        This will be useful for exporting data in non-file formats. For example,
1792        can be used to output grib data directly to S3 using the boto3 client
1793        without the need to write a temporary file to upload first.
1794
1795        Parameters
1796        ----------
1797        validate: default=True
1798            If `True` (DEFAULT), validates first/last four bytes for proper
1799            formatting, else returns None. If `False`, message is output as is.
1800
1801        Returns
1802        -------
1803        to_bytes
1804            Returns GRIB2 formatted message as bytes.
1805        """
1806        if hasattr(self, "_msg"):
1807            if validate:
1808                if self.validate():
1809                    return self._msg
1810                else:
1811                    return None
1812            else:
1813                return self._msg
1814        else:
1815            return None

Return packed GRIB2 message in bytes format.

This will be useful for exporting data in non-file formats. For example, can be used to output grib data directly to S3 using the boto3 client without the need to write a temporary file to upload first.

Parameters
  • validate (default=True): If True (DEFAULT), validates first/last four bytes for proper formatting, else returns None. If False, message is output as is.
Returns
  • to_bytes: Returns GRIB2 formatted message as bytes.
def interpolate( self, method, grid_def_out, method_options=None, drtn=None, num_threads=1):
1817    def interpolate(
1818        self, method, grid_def_out, method_options=None, drtn=None, num_threads=1
1819    ):
1820        """
1821        Grib2Message Interpolator
1822
1823        Performs spatial interpolation via the [NCEPLIBS-ip
1824        library](https://github.com/NOAA-EMC/NCEPLIBS-ip). This interpolate
1825        method only supports scalar interpolation. If you need to perform
1826        vector interpolation, use the module-level `grib2io.interpolate`
1827        function.
1828
1829        Parameters
1830        ----------
1831        method
1832            Interpolate method to use. This can either be an integer or string
1833            using the following mapping:
1834
1835            | Interpolate Scheme | Integer Value |
1836            | :---:              | :---:         |
1837            | 'bilinear'         | 0             |
1838            | 'bicubic'          | 1             |
1839            | 'neighbor'         | 2             |
1840            | 'budget'           | 3             |
1841            | 'spectral'         | 4             |
1842            | 'neighbor-budget'  | 6             |
1843
1844        grid_def_out : grib2io.Grib2GridDef
1845            Grib2GridDef object of the output grid.
1846        method_options : list of ints, optional
1847            Interpolation options. See the NCEPLIBS-ip documentation for
1848            more information on how these are used.
1849        drtn
1850            Data Representation Template to be used for the returned
1851            interpolated GRIB2 message. When `None`, the data representation
1852            template of the source GRIB2 message is used. Once again, it is the
1853            user's responsibility to properly set the Data Representation
1854            Template attributes.
1855        num_threads : int, optional
1856            Number of OpenMP threads to use for interpolation. The default
1857            value is 1. If NCEPLIBS-ip and grib2io's iplib extension module
1858            was not built with OpenMP, then this keyword argument and value
1859            will have no impact.
1860
1861        Returns
1862        -------
1863        interpolate
1864            If interpolating to a grid, a new Grib2Message object is returned.
1865            The GRIB2 metadata of the new Grib2Message object is identical to
1866            the input except where required to be different because of the new
1867            grid specs and possibly a new data representation template.
1868
1869            If interpolating to station points, the interpolated data values are
1870            returned as a numpy.ndarray.
1871        """
1872        section0 = self.section0
1873        section0[-1] = 0
1874        gds = [0, grid_def_out.npoints, 0, 255, grid_def_out.gdtn]
1875        section3 = np.concatenate((gds, grid_def_out.gdt))
1876        drtn = self.drtn if drtn is None else drtn
1877
1878        msg = Grib2Message(
1879            section0,
1880            self.section1,
1881            self.section2,
1882            section3,
1883            self.section4,
1884            None,
1885            self.bitMapFlag.value,
1886            drtn=drtn,
1887        )
1888
1889        msg._msgnum = -1
1890        msg._deflist = self._deflist
1891        msg._coordlist = self._coordlist
1892        if msg.typeOfValues == 0:
1893            pass
1894        elif msg.typeOfValues == 1:
1895            pass
1896        msg._data = interpolate(
1897            self.data,
1898            method,
1899            Grib2GridDef.from_section3(self.section3),
1900            grid_def_out,
1901            method_options=method_options,
1902            num_threads=num_threads,
1903        ).reshape(msg.ny, msg.nx)
1904        msg.section5[0] = grid_def_out.npoints
1905        return msg

Grib2Message Interpolator

Performs spatial interpolation via the NCEPLIBS-ip library. This interpolate method only supports scalar interpolation. If you need to perform vector interpolation, use the module-level grib2io.interpolate function.

Parameters
  • method: Interpolate method to use. This can either be an integer or string using the following mapping:
Interpolate Scheme Integer Value
'bilinear' 0
'bicubic' 1
'neighbor' 2
'budget' 3
'spectral' 4
'neighbor-budget' 6

  • grid_def_out (grib2io.Grib2GridDef): Grib2GridDef object of the output grid.
  • method_options (list of ints, optional): Interpolation options. See the NCEPLIBS-ip documentation for more information on how these are used.
  • drtn: Data Representation Template to be used for the returned interpolated GRIB2 message. When None, the data representation template of the source GRIB2 message is used. Once again, it is the user's responsibility to properly set the Data Representation Template attributes.
  • num_threads (int, optional): Number of OpenMP threads to use for interpolation. The default value is 1. If NCEPLIBS-ip and grib2io's iplib extension module was not built with OpenMP, then this keyword argument and value will have no impact.
Returns
  • interpolate: If interpolating to a grid, a new Grib2Message object is returned. The GRIB2 metadata of the new Grib2Message object is identical to the input except where required to be different because of the new grid specs and possibly a new data representation template.

If interpolating to station points, the interpolated data values are returned as a numpy.ndarray.

def subset(self, lats, lons):
1907    def subset(self, lats, lons):
1908        """
1909        Return a spatial subset.
1910
1911        Currently only supports regular grids of the following types:
1912
1913        | Grid Type                                                    | gdtn  |
1914        | :---:                                                        | :---: |
1915        | Latitude/Longitude, Equidistant Cylindrical, or Plate Carree | 0     |
1916        | Rotated Latitude/Longitude                                   | 1     |
1917        | Mercator                                                     | 10    |
1918        | Polar Stereographic                                          | 20    |
1919        | Lambert Conformal                                            | 30    |
1920        | Albers Equal-Area                                            | 31    |
1921        | Gaussian Latitude/Longitude                                  | 40    |
1922        | Equatorial Azimuthal Equidistant Projection                  | 110   |
1923
1924        Parameters
1925        ----------
1926        lats
1927            List or tuple of latitudes.  The minimum and maximum latitudes will
1928            be used to define the southern and northern boundaries.
1929
1930            The order of the latitudes is not important.  The function will
1931            determine which is the minimum and maximum.
1932
1933            The latitudes should be in decimal degrees with 0.0 at the equator,
1934            positive values in the northern hemisphere increasing to 90, and
1935            negative values in the southern hemisphere decreasing to -90.
1936        lons
1937            List or tuple of longitudes.  The minimum and maximum longitudes
1938            will be used to define the western and eastern boundaries.
1939
1940            The order of the longitudes is not important.  The function will
1941            determine which is the minimum and maximum.
1942
1943            GRIB2 longitudes should be in decimal degrees with 0.0 at the prime
1944            meridian, positive values increasing eastward to 360.  There are no
1945            negative GRIB2 longitudes.
1946
1947            The typical west longitudes that start at 0.0 at the prime meridian
1948            and decrease to -180 westward, are converted to GRIB2 longitudes by
1949            '360 - (absolute value of the west longitude)' where typical
1950            eastern longitudes are unchanged as GRIB2 longitudes.
1951
1952        Returns
1953        -------
1954        subset
1955            A spatial subset of a GRIB2 message.
1956        """
1957        if self.gdtn not in [0, 1, 10, 20, 30, 31, 40, 110]:
1958            raise ValueError(
1959                """
1960
1961Subset only works for
1962    Latitude/Longitude, Equidistant Cylindrical, or Plate Carree (gdtn=0)
1963    Rotated Latitude/Longitude (gdtn=1)
1964    Mercator (gdtn=10)
1965    Polar Stereographic (gdtn=20)
1966    Lambert Conformal (gdtn=30)
1967    Albers Equal-Area (gdtn=31)
1968    Gaussian Latitude/Longitude (gdtn=40)
1969    Equatorial Azimuthal Equidistant Projection (gdtn=110)
1970
1971"""
1972            )
1973
1974        if self.nx == 0 or self.ny == 0:
1975            raise ValueError(
1976                """
1977
1978Subset only works for regular grids.
1979
1980"""
1981            )
1982
1983        newmsg = Grib2Message(
1984            np.copy(self.section0),
1985            np.copy(self.section1),
1986            np.copy(self.section2),
1987            np.copy(self.section3),
1988            np.copy(self.section4),
1989            np.copy(self.section5),
1990        )
1991
1992        msglats, msglons = self.grid()
1993
1994        la1 = np.max(lats)
1995        lo1 = np.min(lons)
1996        la2 = np.min(lats)
1997        lo2 = np.max(lons)
1998
1999        # Find the indices of the first and last grid points to the nearest
2000        # lat/lon values in the grid.
2001        first_lat = np.abs(msglats - la1)
2002        first_lon = np.abs(msglons - lo1)
2003        max_idx = np.maximum(first_lat, first_lon)
2004        first_j, first_i = np.where(max_idx == np.min(max_idx))
2005
2006        last_lat = np.abs(msglats - la2)
2007        last_lon = np.abs(msglons - lo2)
2008        max_idx = np.maximum(last_lat, last_lon)
2009        last_j, last_i = np.where(max_idx == np.min(max_idx))
2010
2011        setattr(newmsg, "latitudeFirstGridpoint", msglats[first_j[0], first_i[0]])
2012        setattr(newmsg, "longitudeFirstGridpoint", msglons[first_j[0], first_i[0]])
2013        setattr(newmsg, "nx", np.abs(first_i[0] - last_i[0]))
2014        setattr(newmsg, "ny", np.abs(first_j[0] - last_j[0]))
2015
2016        # Set *LastGridpoint attributes even if only used for gdtn=[0, 1, 40].
2017        # This information is used to subset xarray datasets and even though
2018        # unnecessary for some supported grid types, it won't affect a grib2io
2019        # message to set them.
2020        setattr(newmsg, "latitudeLastGridpoint", msglats[last_j[0], last_i[0]])
2021        setattr(newmsg, "longitudeLastGridpoint", msglons[last_j[0], last_i[0]])
2022
2023        setattr(
2024            newmsg,
2025            "data",
2026            self.data[
2027                min(first_j[0], last_j[0]) : max(first_j[0], last_j[0]),
2028                min(first_i[0], last_i[0]) : max(first_i[0], last_i[0]),
2029            ].copy(),
2030        )
2031
2032        # Need to set the newmsg._sha1_section3 to a blank string so the grid
2033        # method ignores the cached lat/lon values.  This will force the grid
2034        # method to recompute the lat/lon values for the subsetted grid.
2035        newmsg._sha1_section3 = ""
2036        newmsg.grid()
2037
2038        return newmsg

Return a spatial subset.

Currently only supports regular grids of the following types:

Grid Type gdtn
Latitude/Longitude, Equidistant Cylindrical, or Plate Carree 0
Rotated Latitude/Longitude 1
Mercator 10
Polar Stereographic 20
Lambert Conformal 30
Albers Equal-Area 31
Gaussian Latitude/Longitude 40
Equatorial Azimuthal Equidistant Projection 110
Parameters
  • lats: List or tuple of latitudes. The minimum and maximum latitudes will be used to define the southern and northern boundaries.

The order of the latitudes is not important. The function will determine which is the minimum and maximum.

The latitudes should be in decimal degrees with 0.0 at the equator, positive values in the northern hemisphere increasing to 90, and negative values in the southern hemisphere decreasing to -90.

  • lons: List or tuple of longitudes. The minimum and maximum longitudes will be used to define the western and eastern boundaries.

The order of the longitudes is not important. The function will determine which is the minimum and maximum.

GRIB2 longitudes should be in decimal degrees with 0.0 at the prime meridian, positive values increasing eastward to 360. There are no negative GRIB2 longitudes.

The typical west longitudes that start at 0.0 at the prime meridian and decrease to -180 westward, are converted to GRIB2 longitudes by '360 - (absolute value of the west longitude)' where typical eastern longitudes are unchanged as GRIB2 longitudes.

Returns
  • subset: A spatial subset of a GRIB2 message.
def validate(self):
2040    def validate(self):
2041        """
2042        Validate a complete GRIB2 message.
2043
2044        The g2c library does its own internal validation when g2_gribend() is called, but
2045        we will check in grib2io also. The validation checks if the first 4 bytes in
2046        self._msg is 'GRIB' and '7777' as the last 4 bytes and that the message length in
2047        section 0 equals the length of the packed message.
2048
2049        Returns
2050        -------
2051        `True` if the packed GRIB2 message is complete and well-formed, `False` otherwise.
2052        """
2053        valid = False
2054        if hasattr(self, "_msg"):
2055            if self._msg[0:4] + self._msg[-4:] == b"GRIB7777":
2056                if self.section0[-1] == len(self._msg):
2057                    valid = True
2058        return valid

Validate a complete GRIB2 message.

The g2c library does its own internal validation when g2_gribend() is called, but we will check in grib2io also. The validation checks if the first 4 bytes in self._msg is 'GRIB' and '7777' as the last 4 bytes and that the message length in section 0 equals the length of the packed message.

Returns
  • True if the packed GRIB2 message is complete and well-formed, False otherwise.
@dataclass
class Grib2GridDef:
2590@dataclass
2591class Grib2GridDef:
2592    """
2593    Class for Grid Definition Template Number and Template as attributes.
2594
2595    This allows for cleaner looking code when passing these metadata around.
2596    For example, the `grib2io._Grib2Message.interpolate` method and
2597    `grib2io.interpolate` function accepts these objects.
2598    """
2599
2600    gdtn: int
2601    gdt: NDArray
2602
2603    @classmethod
2604    def from_section3(cls, section3):
2605        return cls(section3[4], section3[5:])
2606
2607    @property
2608    def nx(self):
2609        """Number of grid points in x-direction."""
2610        return int(self.gdt[7])
2611
2612    @property
2613    def ny(self):
2614        """Number of grid points in y-direction."""
2615        return int(self.gdt[8])
2616
2617    @property
2618    def npoints(self):
2619        """Total number of grid points."""
2620        return int(self.gdt[7] * self.gdt[8])
2621
2622    @property
2623    def shape(self):
2624        """Shape of the grid."""
2625        return (int(self.ny), int(self.nx))
2626
2627    def to_section3(self):
2628        """Return a full GRIB2 section3 array."""
2629        return np.array([0, self.npoints, 0, 0, self.gdtn] + list(self.gdt)).astype(
2630            np.int64
2631        )

Class for Grid Definition Template Number and Template as attributes.

This allows for cleaner looking code when passing these metadata around. For example, the grib2io._Grib2Message.interpolate method and grib2io.interpolate function accepts these objects.

Grib2GridDef( gdtn: int, gdt: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]])
gdtn: int
gdt: numpy.ndarray[tuple[typing.Any, ...], numpy.dtype[~_ScalarT]]
@classmethod
def from_section3(cls, section3):
2603    @classmethod
2604    def from_section3(cls, section3):
2605        return cls(section3[4], section3[5:])
nx
2607    @property
2608    def nx(self):
2609        """Number of grid points in x-direction."""
2610        return int(self.gdt[7])

Number of grid points in x-direction.

ny
2612    @property
2613    def ny(self):
2614        """Number of grid points in y-direction."""
2615        return int(self.gdt[8])

Number of grid points in y-direction.

npoints
2617    @property
2618    def npoints(self):
2619        """Total number of grid points."""
2620        return int(self.gdt[7] * self.gdt[8])

Total number of grid points.

shape
2622    @property
2623    def shape(self):
2624        """Shape of the grid."""
2625        return (int(self.ny), int(self.nx))

Shape of the grid.

def to_section3(self):
2627    def to_section3(self):
2628        """Return a full GRIB2 section3 array."""
2629        return np.array([0, self.npoints, 0, 0, self.gdtn] + list(self.gdt)).astype(
2630            np.int64
2631        )

Return a full GRIB2 section3 array.

def msgs_from_index(index: dict, filehandle=None):
836def msgs_from_index(index: dict, filehandle=None):
837    """
838    Construct a list of Grib2Message objects from an index dictionary.
839
840    This function reconstructs a sequence of `Grib2Message` instances using
841    metadata sections stored in an index dictionary. If an open file handle is
842    provided, each message is linked to its on-disk binary data through a
843    `Grib2MessageOnDiskArray`, allowing deferred reading of the actual data
844    values from the GRIB2 file.
845
846    Parameters
847    ----------
848    index : dict
849        Dictionary containing parsed GRIB2 index information, including
850        section data arrays such as ``section0`` through ``section5``,
851        ``sectionOffset``, ``offset``, and ``bmapflag``.
852    filehandle : file-like object, optional
853        An open binary file handle to the GRIB2 file corresponding to the index.
854        If provided, the returned messages can access on-disk data arrays via
855        memory offsets. If not provided, only metadata will be available.
856
857    Returns
858    -------
859    list of Grib2Message
860        List of reconstructed `Grib2Message` objects built from the provided
861        index. Each message contains metadata, and if `filehandle` is given,
862        also references to on-disk data through a `Grib2MessageOnDiskArray`.
863
864    Notes
865    -----
866    - Each message is constructed by zipping the corresponding section entries
867      (sections 0–5 and bitmap flags).
868    - When a file handle is supplied, each message’s `_ondiskarray` attribute is
869      initialized to allow direct access to GRIB2 data values without loading
870      them fully into memory.
871    - The `_msgnum` attribute of each message is assigned sequentially to
872      preserve message order.
873    """
874    zipped = zip(
875        index["section0"],
876        index["section1"],
877        index["section2"],
878        index["section3"],
879        index["section4"],
880        index["section5"],
881        index["bmapflag"],
882    )
883    msgs = [Grib2Message(*sections) for sections in zipped]
884
885    if filehandle is not None:
886        for n, (msg, offset, secpos) in enumerate(
887            zip(msgs, index["offset"], index["sectionOffset"])
888        ):
889            msg._ondiskarray = Grib2MessageOnDiskArray(
890                shape=(msg.ny, msg.nx),
891                ndim=2,
892                dtype=TYPE_OF_VALUES_DTYPE[msg.typeOfValues],
893                filehandle=filehandle,
894                msg=msg,
895                offset=offset,
896                bitmap_offset=secpos[6],
897                data_offset=secpos[7],
898            )
899            msg._msgnum = n
900    return msgs

Construct a list of Grib2Message objects from an index dictionary.

This function reconstructs a sequence of Grib2Message instances using metadata sections stored in an index dictionary. If an open file handle is provided, each message is linked to its on-disk binary data through a Grib2MessageOnDiskArray, allowing deferred reading of the actual data values from the GRIB2 file.

Parameters
  • index (dict): Dictionary containing parsed GRIB2 index information, including section data arrays such as section0 through section5, sectionOffset, offset, and bmapflag.
  • filehandle (file-like object, optional): An open binary file handle to the GRIB2 file corresponding to the index. If provided, the returned messages can access on-disk data arrays via memory offsets. If not provided, only metadata will be available.
Returns
  • list of Grib2Message: List of reconstructed Grib2Message objects built from the provided index. Each message contains metadata, and if filehandle is given, also references to on-disk data through a Grib2MessageOnDiskArray.
Notes
  • Each message is constructed by zipping the corresponding section entries (sections 0–5 and bitmap flags).
  • When a file handle is supplied, each message’s _ondiskarray attribute is initialized to allow direct access to GRIB2 data values without loading them fully into memory.
  • The _msgnum attribute of each message is assigned sequentially to preserve message order.