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

Initialize GRIB2 File object instance.

Parameters
  • filename: File path containing GRIB2 messages OR bytes.
  • 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.

mode
indexfile
save_index
use_index
closed
levels
321    @property
322    def levels(self):
323        """Provides a unique tuple of level strings."""
324        if self._hasindex:
325            return tuple(sorted(set([msg.level for msg in self._msgs])))
326        else:
327            return None

Provides a unique tuple of level strings.

variables
330    @property
331    def variables(self):
332        """Provides a unique tuple of variable shortName strings."""
333        if self._hasindex:
334            return tuple(sorted(set([msg.shortName for msg in self._msgs])))
335        else:
336            return None

Provides a unique tuple of variable shortName strings.

def close(self):
339    def close(self):
340        """Close the file handle."""
341        if not self._filehandle.closed:
342            self.messages = 0
343            self.current_message = 0
344            self._filehandle.close()
345            self.closed = self._filehandle.closed

Close the file handle.

def read(self, size: Optional[int] = None):
348    def read(self, size: Optional[int]=None):
349        """
350        Read size amount of GRIB2 messages from the current position.
351
352        If no argument is given, then size is None and all messages are returned
353        from the current position in the file. This read method follows the
354        behavior of Python's builtin open() function, but whereas that operates
355        on units of bytes, we operate on units of GRIB2 messages.
356
357        Parameters
358        ----------
359        size: default=None
360            The number of GRIB2 messages to read from the current position. If
361            no argument is give, the default value is None and remainder of
362            the file is read.
363
364        Returns
365        -------
366        read
367            ``Grib2Message`` object when size = 1 or a list of Grib2Messages
368            when size > 1.
369        """
370        if size is not None and size < 0:
371            size = None
372        if size is None or size > 1:
373            start = self.tell()
374            stop = self.messages if size is None else start+size
375            if size is None:
376                self.current_message = self.messages-1
377            else:
378                self.current_message += size
379            return self._msgs[slice(start,stop,1)]
380        elif size == 1:
381            self.current_message += 1
382            return self._msgs[self.current_message]
383        else:
384            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):
387    def seek(self, pos: int):
388        """
389        Set the position within the file in units of GRIB2 messages.
390
391        Parameters
392        ----------
393        pos
394            The GRIB2 Message number to set the file pointer to.
395        """
396        if self._hasindex:
397            self._filehandle.seek(self._index['sectionOffset'][0][pos])
398            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):
401    def tell(self):
402        """Returns the position of the file in units of GRIB2 Messages."""
403        return self.current_message

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

def select(self, **kwargs):
406    def select(self, **kwargs):
407        """Select GRIB2 messages by `Grib2Message` attributes."""
408        # TODO: Added ability to process multiple values for each keyword (attribute)
409        idxs = []
410        nkeys = len(kwargs.keys())
411        for k,v in kwargs.items():
412            for m in self._msgs:
413                if hasattr(m,k) and getattr(m,k) == v: idxs.append(m._msgnum)
414        idxs = np.array(idxs,dtype=np.int32)
415        return [self._msgs[i] for i in [ii[0] for ii in collections.Counter(idxs).most_common() if ii[1] == nkeys]]

Select GRIB2 messages by Grib2Message attributes.

def write(self, msg):
418    def write(self, msg):
419        """
420        Writes GRIB2 message object to file.
421
422        Parameters
423        ----------
424        msg
425            GRIB2 message objects to write to file.
426        """
427        if isinstance(msg,list):
428            for m in msg:
429                self.write(m)
430            return
431
432        if issubclass(msg.__class__,_Grib2Message):
433            # TODO: We can consider letting pack return packed bytes instead of associating with message object
434            if hasattr(msg,'_msg'):
435                # write already packed bytes
436                self._filehandle.write(msg._msg)
437            else:
438                if msg._signature == msg._generate_signature() and \
439                   msg._data is None and \
440                   hasattr(msg._ondiskarray,'filehandle'):
441                   # write unchanged message from input
442                   offset = msg._ondiskarray.filehandle.tell()
443                   msg._ondiskarray.filehandle.seek(msg._ondiskarray.offset)
444                   self._filehandle.write(msg._ondiskarray.filehandle.read(msg.section0[-1]))
445                   msg._ondiskarray.filehandle.seek(offset)
446                else:
447                    msg.pack()
448                    self._filehandle.write(msg._msg)
449            self.size = os.path.getsize(self.name)
450            self.messages += 1
451        else:
452            raise TypeError("msg must be a Grib2Message object.")
453        return

Writes GRIB2 message object to file.

Parameters
  • msg: GRIB2 message objects to write to file.
def flush(self):
456    def flush(self):
457        """Flush the file object buffer."""
458        self._filehandle.flush()

Flush the file object buffer.

def levels_by_var(self, name: str):
461    def levels_by_var(self, name: str):
462        """
463        Return a list of level strings given a variable shortName.
464
465        Parameters
466        ----------
467        name
468            Grib2Message variable shortName
469
470        Returns
471        -------
472        levels_by_var
473            A list of unique level strings.
474        """
475        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):
478    def vars_by_level(self, level: str):
479        """
480        Return a list of variable shortName strings given a level.
481
482        Parameters
483        ----------
484        level
485            Grib2Message variable level
486
487        Returns
488        -------
489        vars_by_level
490            A list of unique variable shortName strings.
491        """
492        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
size
def show_config():
34def show_config():
35    """Print grib2io build configuration information."""
36    print(f'grib2io version {__version__} Configuration:')
37    print(f'')
38    print(f'NCEPLIBS-g2c library version: {__g2clib_version__}')
39    print(f'\tStatic library: {g2c_static}')
40    print(f'\tJPEG compression support: {has_jpeg_support}')
41    print(f'\tPNG compression support: {has_png_support}')
42    print(f'\tAEC compression support: {has_aec_support}')
43    print(f'')
44    print(f'NCEPLIPS-ip support: {has_interpolation}')
45    print(f'\tStatic library: {ip_static}')
46    print(f'\tOpenMP support: {has_openmp_support}')
47    print(f'')
48    print(f'Static libs:')
49    for lib in extra_objects:
50        print(f'\t{lib}')
51    print(f'')
52    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):
2035def interpolate(a,
2036    method: Union[int, str],
2037    grid_def_in,
2038    grid_def_out,
2039    method_options=None,
2040    num_threads=1,
2041):
2042    """
2043    This is the module-level interpolation function.
2044
2045    This interfaces with the 4-byte (32-bit float) [NCEPLIBS-ip library](https://github.com/NOAA-EMC/NCEPLIBS-ip)
2046    through grib2io's internal iplib Cython extension module.
2047
2048    Parameters
2049    ----------
2050    a : numpy.ndarray or tuple
2051        Input data.  If `a` is a `numpy.ndarray`, scalar interpolation will be
2052        performed.  If `a` is a `tuple`, then vector interpolation will be
2053        performed with the assumption that u = a[0] and v = a[1] and are both
2054        `numpy.ndarray`.
2055
2056        These data are expected to be in 2-dimensional form with shape (ny, nx)
2057        or 3-dimensional (:, ny, nx) where the 1st dimension represents another
2058        spatial, temporal, or classification (i.e. ensemble members) dimension.
2059        The function will properly flatten the (ny,nx) dimensions into (nx * ny)
2060        acceptable for input into the interpolation subroutines. If needed, these
2061        data will be converted to `np.float32`.
2062    method
2063        Interpolate method to use. This can either be an integer or string using
2064        the following mapping:
2065
2066        | Interpolate Scheme | Integer Value |
2067        | :---:              | :---:         |
2068        | 'bilinear'         | 0             |
2069        | 'bicubic'          | 1             |
2070        | 'neighbor'         | 2             |
2071        | 'budget'           | 3             |
2072        | 'spectral'         | 4             |
2073        | 'neighbor-budget'  | 6             |
2074
2075    grid_def_in : grib2io.Grib2GridDef
2076        Grib2GridDef object for the input grid.
2077    grid_def_out : grib2io.Grib2GridDef
2078        Grib2GridDef object for the output grid or station points.
2079    method_options : list of ints, optional
2080        Interpolation options. See the NCEPLIBS-ip documentation for
2081        more information on how these are used.
2082    num_threads : int, optional
2083        Number of OpenMP threads to use for interpolation. The default
2084        value is 1. If NCEPLIBS-ip and grib2io's iplib extension module
2085        was not built with OpenMP, then this keyword argument and value
2086        will have no impact.
2087
2088    Returns
2089    -------
2090    interpolate
2091        Returns a `numpy.ndarray` of dtype `np.float32` when scalar interpolation
2092        is performed or a `tuple` of `numpy.ndarray`s when vector interpolation is
2093        performed with the assumptions that 0-index is the interpolated u and
2094        1-index is the interpolated v.
2095    """
2096
2097    from . import iplib
2098
2099    prev_num_threads = 1
2100    try:
2101        prev_num_threads = iplib.openmp_get_num_threads()
2102        iplib.openmp_set_num_threads(num_threads)
2103    except(AttributeError):
2104        pass
2105
2106    print(f"grib2io.interpolate thread report: OpenMP num threads = {iplib.openmp_get_num_threads()}")
2107
2108    if isinstance(method,int) and method not in _interp_schemes.values():
2109        raise ValueError('Invalid interpolation method.')
2110    elif isinstance(method,str):
2111        if method in _interp_schemes.keys():
2112            method = _interp_schemes[method]
2113        else:
2114            raise ValueError('Invalid interpolation method.')
2115
2116    if method_options is None:
2117        method_options = np.zeros((20),dtype=np.int32)
2118        if method in {3,6}:
2119            method_options[0:2] = -1
2120
2121    mi = grid_def_in.npoints
2122    mo = grid_def_out.npoints
2123
2124    # Adjust shape of input array(s)
2125    a, newshp = _adjust_array_shape_for_interp(a, grid_def_in, grid_def_out)
2126
2127    # Call interpolation subroutines according to type of a.
2128    if isinstance(a,np.ndarray):
2129        # Scalar
2130        km = a.shape[0]
2131        if np.any(np.isnan(a)):
2132            ibi = np.ones((km), dtype=np.int32)
2133            li = np.where(np.isnan(a),0,1).astype(np.uint8)
2134        else:
2135            ibi = np.zeros((km), dtype=np.int32)
2136            li = np.zeros(a.shape,dtype=np.uint8)
2137        no, rlat, rlon, ibo, lo, go, iret = iplib.interpolate_scalar(
2138            method,
2139            method_options,
2140            grid_def_in.gdtn,
2141            np.array(grid_def_in.gdt, dtype=np.int32),
2142            grid_def_out.gdtn,
2143            np.array(grid_def_out.gdt, dtype=np.int32),
2144            mi,
2145            mo,
2146            km,
2147            ibi,
2148            li,
2149            a,
2150        )
2151        out = np.where(lo==0,np.nan,go).reshape(newshp)
2152    elif isinstance(a,tuple):
2153        # Vector
2154        km = a[0].shape[0]
2155        if np.any(np.isnan(a)):
2156            ibi = np.ones((km), dtype=np.int32)
2157            li = np.where(np.isnan(a),0,1).astype(np.uint8)
2158        else:
2159            ibi = np.zeros((km), dtype=np.int32)
2160            li = np.zeros(a[0].shape,dtype=np.uint8)
2161        no, rlat, rlon, crot, srot, ibo, lo, uo, vo, iret = iplib.interpolate_vector(
2162            method,
2163            method_options,
2164            grid_def_in.gdtn,
2165            np.array(grid_def_in.gdt, dtype=np.int32),
2166            grid_def_out.gdtn,
2167            np.array(grid_def_out.gdt, dtype=np.int32),
2168            mi,
2169            mo,
2170            km,
2171            ibi,
2172            li,
2173            a[0],
2174            a[1],
2175        )
2176        uo = np.where(lo==0,np.nan,uo).reshape(newshp)
2177        vo = np.where(lo==0,np.nan,vo).reshape(newshp)
2178        out = (uo,vo)
2179
2180    try:
2181        iplib.openmp_set_num_threads(prev_num_threads)
2182    except(AttributeError):
2183        pass
2184
2185    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):
2188def interpolate_to_stations(
2189    a,
2190    method: Union[int, str],
2191    grid_def_in,
2192    lats,
2193    lons,
2194    method_options=None,
2195    num_threads=1
2196):
2197    """
2198    Module-level interpolation function for interpolation to stations.
2199
2200    Interfaces with the 4-byte (32-bit float) [NCEPLIBS-ip library](https://github.com/NOAA-EMC/NCEPLIBS-ip)
2201    via grib2io's iplib Cython exntension module. It supports scalar and
2202    vector interpolation according to the type of object a.
2203
2204    Parameters
2205    ----------
2206    a : numpy.ndarray or tuple
2207        Input data.  If `a` is a `numpy.ndarray`, scalar interpolation will be
2208        performed.  If `a` is a `tuple`, then vector interpolation will be
2209        performed with the assumption that u = a[0] and v = a[1] and are both
2210        `numpy.ndarray`.
2211
2212        These data are expected to be in 2-dimensional form with shape (ny, nx)
2213        or 3-dimensional (:, ny, nx) where the 1st dimension represents another
2214        spatial, temporal, or classification (i.e. ensemble members) dimension.
2215        The function will properly flatten the (ny,nx) dimensions into (nx * ny)
2216        acceptable for input into the interpolation subroutines. If needed, these
2217        data will be converted to `np.float32`.
2218    method
2219        Interpolate method to use. This can either be an integer or string using
2220        the following mapping:
2221
2222        | Interpolate Scheme | Integer Value |
2223        | :---:              | :---:         |
2224        | 'bilinear'         | 0             |
2225        | 'bicubic'          | 1             |
2226        | 'neighbor'         | 2             |
2227        | 'budget'           | 3             |
2228        | 'spectral'         | 4             |
2229        | 'neighbor-budget'  | 6             |
2230
2231    grid_def_in : grib2io.Grib2GridDef
2232        Grib2GridDef object for the input grid.
2233    lats : numpy.ndarray or list
2234        Latitudes for station points
2235    lons : numpy.ndarray or list
2236        Longitudes for station points
2237    method_options : list of ints, optional
2238        Interpolation options. See the NCEPLIBS-ip documentation for
2239        more information on how these are used.
2240    num_threads : int, optional
2241        Number of OpenMP threads to use for interpolation. The default
2242        value is 1. If NCEPLIBS-ip and grib2io's iplib extension module
2243        was not built with OpenMP, then this keyword argument and value
2244        will have no impact.
2245
2246    Returns
2247    -------
2248    interpolate_to_stations
2249        Returns a `numpy.ndarray` of dtype `np.float32` when scalar
2250        interpolation is performed or a `tuple` of `numpy.ndarray`s
2251        when vector interpolation is performed with the assumptions
2252        that 0-index is the interpolated u and 1-index is the
2253        interpolated v.
2254    """
2255    from . import iplib
2256
2257    # Define function to apply mask when stations are outside grid domain
2258    def _reshape_and_mask_post_interp(a, shape, mask):
2259        a = a.reshape(shape)
2260        if a.shape[-1] != mask.shape[0]:
2261            raise ValueError(f"Station mask length does not match interpolated data.")
2262        a[..., mask] = np.nan
2263        return a
2264
2265    prev_num_threads = 1
2266    try:
2267        prev_num_threads = iplib.openmp_get_num_threads()
2268        iplib.openmp_set_num_threads(num_threads)
2269    except(AttributeError):
2270        pass
2271
2272    if isinstance(method,int) and method not in _interp_schemes.values():
2273        raise ValueError('Invalid interpolation method.')
2274    elif isinstance(method,str):
2275        if method in _interp_schemes.keys():
2276            method = _interp_schemes[method]
2277        else:
2278            raise ValueError('Invalid interpolation method.')
2279
2280    if method_options is None:
2281        method_options = np.zeros((20),dtype=np.int32)
2282        if method in {3,6}:
2283            method_options[0:2] = -1
2284
2285    # Check lats and lons
2286    if isinstance(lats,list):
2287        nlats = len(lats)
2288    elif isinstance(lats,np.ndarray) and len(lats.shape) == 1:
2289        nlats = lats.shape[0]
2290    else:
2291        raise ValueError('Station latitudes must be a list or 1-D NumPy array.')
2292    if isinstance(lons,list):
2293        nlons = len(lons)
2294    elif isinstance(lons,np.ndarray) and len(lons.shape) == 1:
2295        nlons = lons.shape[0]
2296    else:
2297        raise ValueError('Station longitudes must be a list or 1-D NumPy array.')
2298    if nlats != nlons:
2299        raise ValueError('Station lats and lons must be same size.')
2300
2301    mi = grid_def_in.npoints
2302    mo = nlats
2303
2304    # Adjust shape of input array(s)
2305    a, newshp = _adjust_array_shape_for_interp_stations(a, grid_def_in, mo)
2306
2307    # Use gdtn = -1 for stations and an empty template array
2308    gdtn_out = -1
2309    gdt_out = np.zeros((200),dtype=np.int32)
2310
2311    # Before we interpolate, get the grid coordinates for stations.
2312    xloc, yloc = utils.latlon_to_ij(
2313        grid_def_in.gdtn,
2314        grid_def_in.gdt,
2315        np.array(lats, dtype=np.float32),
2316        np.array(lons, dtype=np.float32),
2317    )
2318    xloc_mask = np.where(np.isnan(xloc), True, False)
2319    yloc_mask = np.where(np.isnan(yloc), True, False)
2320    mask = xloc_mask & yloc_mask
2321
2322    # Call interpolation subroutines according to type of a.
2323    if isinstance(a,np.ndarray):
2324        # Scalar
2325        km = a.shape[0]
2326        ibi = np.zeros((km), dtype=np.int32)
2327        li = np.zeros(a.shape,dtype=np.uint8)
2328        no, rlat, rlon, ibo, lo, go, iret = iplib.interpolate_scalar(
2329            method,
2330            method_options,
2331            grid_def_in.gdtn,
2332            np.array(grid_def_in.gdt, dtype=np.int32),
2333            gdtn_out,
2334            gdt_out,
2335            mi,
2336            mo,
2337            km,
2338            ibi,
2339            li,
2340            a,
2341            lats=np.array(lats, dtype=np.float32),
2342            lons=np.array(lons, dtype=np.float32),
2343        )
2344        out = _reshape_and_mask_post_interp(go, newshp, mask)
2345
2346    elif isinstance(a,tuple):
2347        # Vector
2348        km = a[0].shape[0]
2349        ibi = np.zeros((km), dtype=np.int32)
2350        li = np.zeros(a[0].shape,dtype=np.uint8)
2351        no, rlat, rlon, crot, srot, ibo, lo, uo, vo, iret = iplib.interpolate_vector(
2352            method,
2353            method_options,
2354            grid_def_in.gdtn,
2355            np.array(grid_def_in.gdt, dtype=np.int32),
2356            gdtn_out,
2357            gdt_out,
2358            mi,
2359            mo,
2360            km,
2361            ibi,
2362            li,
2363            a[0],
2364            a[1],
2365            lats=np.array(lats, dtype=np.float32),
2366            lons=np.array(lons, dtype=np.float32),
2367        )
2368        out = (
2369            _reshape_and_mask_post_interp(uo, newshp, mask),
2370            _reshape_and_mask_post_interp(vo, newshp, mask)
2371        )
2372
2373    try:
2374        iplib.openmp_set_num_threads(prev_num_threads)
2375    except(AttributeError):
2376        pass
2377
2378    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:
839class Grib2Message:
840    """
841    Creation class for a GRIB2 message.
842
843    This class returns a dynamically-created Grib2Message object that
844    inherits from `_Grib2Message` and grid, product, data representation
845    template classes according to the template numbers for the respective
846    sections. If `section3`, `section4`, or `section5` are omitted, then
847    the appropriate keyword arguments for the template number `gdtn=`,
848    `pdtn=`, or `drtn=` must be provided.
849
850    Parameters
851    ----------
852    section0
853        GRIB2 section 0 array.
854    section1
855        GRIB2 section 1 array.
856    section2
857        Local Use section data.
858    section3
859        GRIB2 section 3 array.
860    section4
861        GRIB2 section 4 array.
862    section5
863        GRIB2 section 5 array.
864
865    Returns
866    -------
867    Msg
868        A dynamically-create Grib2Message object that inherits from
869        _Grib2Message, a grid definition template class, product
870        definition template class, and a data representation template
871        class.
872    """
873    def __new__(self, section0: NDArray = np.array([struct.unpack('>I',b'GRIB')[0],0,0,2,0]),
874                      section1: NDArray = np.zeros((13),dtype=np.int64),
875                      section2: Optional[bytes] = None,
876                      section3: Optional[NDArray] = None,
877                      section4: Optional[NDArray] = None,
878                      section5: Optional[NDArray] = None, *args, **kwargs):
879
880        if np.all(section1==0):
881            try:
882                # Python >= 3.10
883                section1[5:11] = datetime.datetime.fromtimestamp(0, datetime.UTC).timetuple()[:6]
884            except(AttributeError):
885                # Python < 3.10
886                section1[5:11] = datetime.datetime.utcfromtimestamp(0).timetuple()[:6]
887
888        bases = list()
889        if section3 is None:
890            if 'gdtn' in kwargs.keys():
891                gdtn = kwargs['gdtn']
892                Gdt = templates.gdt_class_by_gdtn(gdtn)
893                bases.append(Gdt)
894                section3 = np.zeros((Gdt._len+5),dtype=np.int64)
895                section3[4] = gdtn
896            else:
897                raise ValueError("Must provide GRIB2 Grid Definition Template Number or section 3 array")
898        else:
899            gdtn = section3[4]
900            Gdt = templates.gdt_class_by_gdtn(gdtn)
901            bases.append(Gdt)
902
903        if section4 is None:
904            if 'pdtn' in kwargs.keys():
905                pdtn = kwargs['pdtn']
906                Pdt = templates.pdt_class_by_pdtn(pdtn)
907                bases.append(Pdt)
908                section4 = np.zeros((Pdt._len+2),dtype=np.int64)
909                section4[1] = pdtn
910            else:
911                raise ValueError("Must provide GRIB2 Production Definition Template Number or section 4 array")
912        else:
913            pdtn = section4[1]
914            Pdt = templates.pdt_class_by_pdtn(pdtn)
915            bases.append(Pdt)
916
917        if section5 is None:
918            if 'drtn' in kwargs.keys():
919                drtn = kwargs['drtn']
920                Drt = templates.drt_class_by_drtn(drtn)
921                bases.append(Drt)
922                section5 = np.zeros((Drt._len+2),dtype=np.int64)
923                section5[1] = drtn
924            else:
925                raise ValueError("Must provide GRIB2 Data Representation Template Number or section 5 array")
926        else:
927            drtn = section5[1]
928            Drt = templates.drt_class_by_drtn(drtn)
929            bases.append(Drt)
930
931        # attempt to use existing Msg class if it has already been made with gdtn,pdtn,drtn combo
932        try:
933            Msg = _msg_class_store[f"{gdtn}:{pdtn}:{drtn}"]
934        except KeyError:
935            @dataclass(init=False, repr=False)
936            class Msg(_Grib2Message, *bases):
937                pass
938            _msg_class_store[f"{gdtn}:{pdtn}:{drtn}"] = Msg
939
940        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:
 943@dataclass
 944class _Grib2Message:
 945    """
 946    GRIB2 Message base class.
 947    """
 948    # GRIB2 Sections
 949    section0: NDArray = field(init=True,repr=False)
 950    section1: NDArray = field(init=True,repr=False)
 951    section2: bytes = field(init=True,repr=False)
 952    section3: NDArray = field(init=True,repr=False)
 953    section4: NDArray = field(init=True,repr=False)
 954    section5: NDArray = field(init=True,repr=False)
 955    bitMapFlag: templates.Grib2Metadata = field(init=True,repr=False,default=255)
 956
 957    # Section 0 looked up attributes
 958    indicatorSection: NDArray = field(init=False,repr=False,default=templates.IndicatorSection())
 959    discipline: templates.Grib2Metadata = field(init=False,repr=False,default=templates.Discipline())
 960
 961    # Section 1 looked up attributes
 962    identificationSection: NDArray = field(init=False,repr=False,default=templates.IdentificationSection())
 963    originatingCenter: templates.Grib2Metadata = field(init=False,repr=False,default=templates.OriginatingCenter())
 964    originatingSubCenter: templates.Grib2Metadata = field(init=False,repr=False,default=templates.OriginatingSubCenter())
 965    masterTableInfo: templates.Grib2Metadata = field(init=False,repr=False,default=templates.MasterTableInfo())
 966    localTableInfo: templates.Grib2Metadata = field(init=False,repr=False,default=templates.LocalTableInfo())
 967    significanceOfReferenceTime: templates.Grib2Metadata = field(init=False,repr=False,default=templates.SignificanceOfReferenceTime())
 968    year: int = field(init=False,repr=False,default=templates.Year())
 969    month: int = field(init=False,repr=False,default=templates.Month())
 970    day: int = field(init=False,repr=False,default=templates.Day())
 971    hour: int = field(init=False,repr=False,default=templates.Hour())
 972    minute: int = field(init=False,repr=False,default=templates.Minute())
 973    second: int = field(init=False,repr=False,default=templates.Second())
 974    refDate: datetime.datetime = field(init=False,repr=False,default=templates.RefDate())
 975    productionStatus: templates.Grib2Metadata = field(init=False,repr=False,default=templates.ProductionStatus())
 976    typeOfData: templates.Grib2Metadata = field(init=False,repr=False,default=templates.TypeOfData())
 977
 978    # Section 3 looked up common attributes.  Other looked up attributes are available according
 979    # to the Grid Definition Template.
 980    gridDefinitionSection: NDArray = field(init=False,repr=False,default=templates.GridDefinitionSection())
 981    sourceOfGridDefinition: int = field(init=False,repr=False,default=templates.SourceOfGridDefinition())
 982    numberOfDataPoints: int = field(init=False,repr=False,default=templates.NumberOfDataPoints())
 983    interpretationOfListOfNumbers: templates.Grib2Metadata = field(init=False,repr=False,default=templates.InterpretationOfListOfNumbers())
 984    gridDefinitionTemplateNumber: templates.Grib2Metadata = field(init=False,repr=False,default=templates.GridDefinitionTemplateNumber())
 985    gridDefinitionTemplate: list = field(init=False,repr=False,default=templates.GridDefinitionTemplate())
 986    _earthparams: dict = field(init=False,repr=False,default=templates.EarthParams())
 987    _dxsign: float = field(init=False,repr=False,default=templates.DxSign())
 988    _dysign: float = field(init=False,repr=False,default=templates.DySign())
 989    _llscalefactor: float = field(init=False,repr=False,default=templates.LLScaleFactor())
 990    _lldivisor: float = field(init=False,repr=False,default=templates.LLDivisor())
 991    _xydivisor: float = field(init=False,repr=False,default=templates.XYDivisor())
 992    shapeOfEarth: templates.Grib2Metadata = field(init=False,repr=False,default=templates.ShapeOfEarth())
 993    earthShape: str = field(init=False,repr=False,default=templates.EarthShape())
 994    earthRadius: float = field(init=False,repr=False,default=templates.EarthRadius())
 995    earthMajorAxis: float = field(init=False,repr=False,default=templates.EarthMajorAxis())
 996    earthMinorAxis: float = field(init=False,repr=False,default=templates.EarthMinorAxis())
 997    resolutionAndComponentFlags: list = field(init=False,repr=False,default=templates.ResolutionAndComponentFlags())
 998    ny: int = field(init=False,repr=False,default=templates.Ny())
 999    nx: int = field(init=False,repr=False,default=templates.Nx())
1000    scanModeFlags: list = field(init=False,repr=False,default=templates.ScanModeFlags())
1001    projParameters: dict = field(init=False,repr=False,default=templates.ProjParameters())
1002
1003    # Section 4
1004    productDefinitionTemplateNumber: templates.Grib2Metadata = field(init=False,repr=False,default=templates.ProductDefinitionTemplateNumber())
1005    productDefinitionTemplate: NDArray = field(init=False,repr=False,default=templates.ProductDefinitionTemplate())
1006
1007    # Section 5 looked up common attributes.  Other looked up attributes are
1008    # available according to the Data Representation Template.
1009    numberOfPackedValues: int = field(init=False,repr=False,default=templates.NumberOfPackedValues())
1010    dataRepresentationTemplateNumber: templates.Grib2Metadata = field(init=False,repr=False,default=templates.DataRepresentationTemplateNumber())
1011    dataRepresentationTemplate: list = field(init=False,repr=False,default=templates.DataRepresentationTemplate())
1012    typeOfValues: templates.Grib2Metadata = field(init=False,repr=False,default=templates.TypeOfValues())
1013
1014    def __copy__(self):
1015        """Shallow copy"""
1016        new = Grib2Message(
1017            self.section0,
1018            self.section1,
1019            self.section2,
1020            self.section3,
1021            self.section4,
1022            drtn=self.drtn,
1023        )
1024        return new
1025
1026    def __deepcopy__(self, memo):
1027        """Deep copy"""
1028        new = Grib2Message(
1029            np.copy(self.section0),
1030            np.copy(self.section1),
1031            copy.deepcopy(self.section2),
1032            np.copy(self.section3),
1033            np.copy(self.section4),
1034            np.copy(self.section5),
1035        )
1036        memo[id(self)] = new
1037        new.data = np.copy(self.data)
1038        new.bitmap = None if self.bitmap is None else np.copy(self.bitmap)
1039        return new
1040
1041    def __post_init__(self):
1042        """Set some attributes after init."""
1043        self._auto_nans = _AUTO_NANS
1044        self._coordlist = np.zeros((0), dtype=np.float32)
1045        self._data = None
1046        self._deflist = np.zeros((0), dtype=np.int64)
1047        self._msgnum = -1
1048        self._ondiskarray = None
1049        self._orig_section5 = np.copy(self.section5)
1050        self._signature = self._generate_signature()
1051        try:
1052            self._sha1_section3 = hashlib.sha1(self.section3).hexdigest()
1053        except(TypeError):
1054            pass
1055        self.bitMapFlag = templates.Grib2Metadata(self.bitMapFlag,table='6.0')
1056        self.bitmap = None
1057
1058    @property
1059    def _isNDFD(self):
1060        """Check if GRIB2 message is from NWS NDFD"""
1061        return np.all(self.section1[0:2]==[8,65535])
1062
1063    @property
1064    def _isAerosol(self):
1065        """Check if GRIB2 message contains aerosol data"""
1066        is_aero_template = self.productDefinitionTemplateNumber.value in tables.AEROSOL_PDTNS
1067        is_aero_param = ((str(self.parameterCategory) == '13') |
1068                     (str(self.parameterCategory) == '20')) and str(self.parameterNumber) in tables.AEROSOL_PARAMS
1069        # Check table 4.205 aerosol presence
1070        is_aero_type = (str(self.parameterCategory) == '205' and
1071                        str(self.parameterNumber) == '1')
1072        return is_aero_template or is_aero_param or is_aero_type
1073
1074    @property
1075    def gdtn(self):
1076        """Return Grid Definition Template Number"""
1077        return self.section3[4]
1078
1079
1080    @property
1081    def gdt(self):
1082        """Return Grid Definition Template."""
1083        return self.gridDefinitionTemplate
1084
1085
1086    @property
1087    def pdtn(self):
1088        """Return Product Definition Template Number."""
1089        return self.section4[1]
1090
1091
1092    @property
1093    def pdt(self):
1094        """Return Product Definition Template."""
1095        return self.productDefinitionTemplate
1096
1097
1098    @property
1099    def drtn(self):
1100        """Return Data Representation Template Number."""
1101        return self.section5[1]
1102
1103
1104    @property
1105    def drt(self):
1106        """Return Data Representation Template."""
1107        return self.dataRepresentationTemplate
1108
1109
1110    @property
1111    def pdy(self):
1112        """Return the PDY ('YYYYMMDD')."""
1113        return ''.join([str(i) for i in self.section1[5:8]])
1114
1115
1116    @property
1117    def griddef(self):
1118        """Return a Grib2GridDef instance for a GRIB2 message."""
1119        return Grib2GridDef.from_section3(self.section3)
1120
1121
1122    @property
1123    def lats(self):
1124        """Return grid latitudes."""
1125        return self.latlons()[0]
1126
1127
1128    @property
1129    def lons(self):
1130        """Return grid longitudes."""
1131        return self.latlons()[1]
1132
1133    @property
1134    def min(self):
1135        """Return minimum value of data."""
1136        return np.nanmin(self.data)
1137
1138
1139    @property
1140    def max(self):
1141        """Return maximum value of data."""
1142        return np.nanmax(self.data)
1143
1144
1145    @property
1146    def mean(self):
1147        """Return mean value of data."""
1148        return np.nanmean(self.data)
1149
1150
1151    @property
1152    def median(self):
1153        """Return median value of data."""
1154        return np.nanmedian(self.data)
1155
1156
1157    @property
1158    def shape(self):
1159        """Return shape of data."""
1160        return self.griddef.shape
1161
1162
1163    def __repr__(self):
1164        """
1165        Return an unambiguous string representation of the object.
1166
1167        Returns
1168        -------
1169        repr
1170            A string representation of the object, including information from
1171            sections 0, 1, 3, 4, 5, and 6.
1172        """
1173        info = ''
1174        for sect in [0,1,3,4,5,6]:
1175            for k,v in self.attrs_by_section(sect,values=True).items():
1176                info += f'Section {sect}: {k} = {v}\n'
1177        return info
1178
1179    def __str__(self):
1180        """
1181        Return a readable string representation of the object.
1182
1183        Returns
1184        -------
1185        str
1186            A formatted string representation of the object, including
1187            selected attributes.
1188        """
1189        return (f'{self._msgnum}:d={self.refDate}:{self.shortName}:'
1190                f'{self.fullName} ({self.units}):{self.level}:'
1191                f'{self.leadTime}')
1192
1193
1194    def _generate_signature(self):
1195        """Generature SHA-1 hash string from GRIB2 integer sections."""
1196        return hashlib.sha1(np.concatenate((self.section0,self.section1,
1197                                            self.section3,self.section4,
1198                                            self.section5))).hexdigest()
1199
1200
1201    def attrs_by_section(self, sect: int, values: bool=False):
1202        """
1203        Provide a tuple of attribute names for the given GRIB2 section.
1204
1205        Parameters
1206        ----------
1207        sect
1208            The GRIB2 section number.
1209        values
1210            Optional (default is `False`) argument to return attributes values.
1211
1212        Returns
1213        -------
1214        attrs_by_section
1215            A list of attribute names or dict of name:value pairs if `values =
1216            True`.
1217        """
1218        if sect in {0,1,6}:
1219            attrs = templates._section_attrs[sect]
1220        elif sect in {3,4,5}:
1221            def _find_class_index(n):
1222                _key = {3:'Grid', 4:'Product', 5:'Data'}
1223                for i,c in enumerate(self.__class__.__mro__):
1224                    if _key[n] in c.__name__:
1225                        return i
1226                else:
1227                    return []
1228            if sys.version_info.minor <= 8:
1229                attrs = templates._section_attrs[sect]+\
1230                        [a for a in dir(self.__class__.__mro__[_find_class_index(sect)]) if not a.startswith('_')]
1231            else:
1232                attrs = templates._section_attrs[sect]+\
1233                        self.__class__.__mro__[_find_class_index(sect)]._attrs()
1234        else:
1235            attrs = []
1236        if values:
1237            return {k:getattr(self,k) for k in attrs}
1238        else:
1239            return attrs
1240
1241
1242    def copy(self, deep: bool = True):
1243        """Returns a copy of this Grib2Message.
1244
1245        When `deep=True`, a copy is made of each of the GRIB2 section arrays and
1246        the data are unpacked from the source object and copied into the new
1247        object. Otherwise, a shallow copy of each array is performed and no data
1248        are copied.
1249
1250        Parameters
1251        ----------
1252        deep : bool, default: True
1253            Whether each GRIB2 section array and data are copied onto
1254            the new object. Default is True.
1255
1256        Returns
1257        -------
1258        object : Grib2Message
1259            New Grib2Message object.
1260
1261            .. versionadded:: 2.6.0
1262        """
1263        return copy.deepcopy(self) if deep else copy.copy(self)
1264
1265
1266    def pack(self):
1267        """
1268        Pack GRIB2 section data into a binary message.
1269
1270        It is the user's responsibility to populate the GRIB2 section
1271        information with appropriate metadata.
1272        """
1273        # Create beginning of packed binary message with section 0 and 1 data.
1274        self._sections = []
1275        self._msg,self._pos = g2clib.grib2_create(self.indicatorSection[2:4],self.identificationSection)
1276        self._sections += [0,1]
1277
1278        # Add section 2 if present.
1279        if isinstance(self.section2,bytes) and len(self.section2) > 0:
1280            self._msg,self._pos = g2clib.grib2_addlocal(self._msg,self.section2)
1281            self._sections.append(2)
1282
1283        # Add section 3.
1284        self.section3[1] = self.nx * self.ny
1285        self._msg,self._pos = g2clib.grib2_addgrid(self._msg,self.gridDefinitionSection,
1286                                                   self.gridDefinitionTemplate,
1287                                                   self._deflist)
1288        self._sections.append(3)
1289
1290        # Prepare data.
1291        if self._data is None:
1292            if self._ondiskarray is None:
1293                raise ValueError("Grib2Message object has no data, thus it cannot be packed.")
1294        field = np.copy(self.data)
1295        if self.scanModeFlags is not None:
1296            if self.scanModeFlags[3]:
1297                fieldsave = field.astype('f') # Casting makes a copy
1298                field[1::2,:] = fieldsave[1::2,::-1]
1299        fld = field.astype('f')
1300
1301        # Prepare bitmap, if necessary
1302        bitmapflag = self.bitMapFlag.value
1303        if bitmapflag == 0:
1304            if self.bitmap is not None:
1305                bmap = np.ravel(self.bitmap).astype(DEFAULT_NUMPY_INT)
1306            else:
1307                bmap = np.ravel(np.where(np.isnan(fld),0,1)).astype(DEFAULT_NUMPY_INT)
1308        else:
1309            bmap = None
1310
1311        # Prepare data for packing if nans are present
1312        fld = np.ravel(fld)
1313        if bitmapflag in {0,254}:
1314            fld = np.where(np.isnan(fld),0,fld)
1315        else:
1316            if np.isnan(fld).any():
1317                if hasattr(self,'priMissingValue'):
1318                    fld = np.where(np.isnan(fld),self.priMissingValue,fld)
1319            if hasattr(self,'_missvalmap'):
1320                if hasattr(self,'priMissingValue'):
1321                    fld = np.where(self._missvalmap==1,self.priMissingValue,fld)
1322                if hasattr(self,'secMissingValue'):
1323                    fld = np.where(self._missvalmap==2,self.secMissingValue,fld)
1324
1325        # Add sections 4, 5, 6, and 7.
1326        self._msg,self._pos = g2clib.grib2_addfield(self._msg,self.pdtn,
1327                                                    self.productDefinitionTemplate,
1328                                                    self._coordlist,
1329                                                    self.drtn,
1330                                                    self.dataRepresentationTemplate,
1331                                                    fld,
1332                                                    bitmapflag,
1333                                                    bmap)
1334        self._sections.append(4)
1335        self._sections.append(5)
1336        self._sections.append(6)
1337        self._sections.append(7)
1338
1339        # Finalize GRIB2 message with section 8.
1340        self._msg, self._pos = g2clib.grib2_end(self._msg)
1341        self._sections.append(8)
1342        self.section0[-1] = len(self._msg)
1343
1344
1345    @property
1346    def data(self) -> np.array:
1347        """Access the unpacked data values."""
1348        if self._data is None:
1349            if self._auto_nans != _AUTO_NANS:
1350                self._data = self._ondiskarray
1351            self._data = np.asarray(self._ondiskarray)
1352        return self._data
1353
1354
1355    @data.setter
1356    def data(self, arr):
1357        """
1358        Set the internal data array, enforcing shape (ny, nx) and dtype float32.
1359
1360        If the Grid Definition Section (section 3) of Grib2Message object is
1361        not fully formed (i.e. nx, ny = 0, 0), then the shape of the data array
1362        will be used to set nx and ny of the Grib2Message object. It will be the
1363        responsibility of the user to populate the rest of the Grid Definition
1364        Section attributes.
1365
1366        Parameters
1367        ----------
1368        arr : array_like
1369            A 2D array whose shape must match ``(self.ny, self.nx)``.
1370            It will be converted to ``float32`` and C-contiguous if needed.
1371
1372        Raises
1373        ------
1374        ValueError
1375            If the shape of `arr` does not match the expected dimensions.
1376        """
1377        if not isinstance(arr, np.ndarray):
1378            raise ValueError(f"Grib2Message data only supports numpy arrays")
1379        if self.nx == 0 and self.ny == 0:
1380            self.ny = arr.shape[0]
1381            self.nx = arr.shape[1]
1382        if arr.shape != (self.ny, self.nx):
1383            raise ValueError(
1384                f"Data shape mismatch: expected ({self.ny}, {self.nx}), "
1385                f"got {arr.shape}"
1386            )
1387        # Ensure contiguous memory layout (important for C interoperability)
1388        if not arr.flags["C_CONTIGUOUS"]:
1389            arr = np.ascontiguousarray(arr, dtype=np.float32)
1390        self._data = arr
1391
1392
1393    def flush_data(self):
1394        """
1395        Flush the unpacked data values from the Grib2Message object.
1396
1397        Notes
1398        -----
1399        If the Grib2Message object was constructed from "scratch" (i.e.
1400        not read from file), this method will remove the data array from
1401        the object and it cannot be recovered.
1402        """
1403        self._data = None
1404        self.bitmap = None
1405
1406
1407    def __getitem__(self, item):
1408        return self.data[item]
1409
1410
1411    def __setitem__(self, item):
1412        raise NotImplementedError('assignment of data not supported via setitem')
1413
1414
1415    def latlons(self, *args, **kwrgs):
1416        """Alias for `grib2io.Grib2Message.grid` method."""
1417        return self.grid(*args, **kwrgs)
1418
1419
1420    def grid(self, unrotate: bool=True):
1421        """
1422        Return lats,lons (in degrees) of grid.
1423
1424        Currently can handle reg. lat/lon,cglobal Gaussian, mercator,
1425        stereographic, lambert conformal, albers equal-area, space-view and
1426        azimuthal equidistant grids.
1427
1428        Parameters
1429        ----------
1430        unrotate
1431            If `True` [DEFAULT], and grid is rotated lat/lon, then unrotate the
1432            grid, otherwise `False`, do not.
1433
1434        Returns
1435        -------
1436        lats, lons : numpy.ndarray
1437            Returns two numpy.ndarrays with dtype=numpy.float32 of grid
1438            latitudes and longitudes in units of degrees.
1439        """
1440        if self._sha1_section3 in _latlon_datastore.keys():
1441            return (_latlon_datastore[self._sha1_section3]['latitude'],
1442                    _latlon_datastore[self._sha1_section3]['longitude'])
1443        gdtn = self.gridDefinitionTemplateNumber.value
1444        gdtmpl = self.gridDefinitionTemplate
1445        reggrid = self.gridDefinitionSection[2] == 0 # This means regular 2-d grid
1446        if gdtn == 0:
1447            # Regular lat/lon grid
1448            lon1, lat1 = self.longitudeFirstGridpoint, self.latitudeFirstGridpoint
1449            lon2, lat2 = self.longitudeLastGridpoint, self.latitudeLastGridpoint
1450            dlon = self.gridlengthXDirection
1451            dlat = self.gridlengthYDirection
1452            if lon2 < lon1 and dlon < 0: lon1 = -lon1
1453            lats = np.linspace(lat1,lat2,self.ny)
1454            if reggrid:
1455                lons = np.linspace(lon1,lon2,self.nx)
1456            else:
1457                lons = np.linspace(lon1,lon2,self.ny*2)
1458            lons,lats = np.meshgrid(lons,lats) # Make 2-d arrays.
1459        elif gdtn == 1: # Rotated Lat/Lon grid
1460            pj = pyproj.Proj(self.projParameters)
1461            lat1,lon1 = self.latitudeFirstGridpoint,self.longitudeFirstGridpoint
1462            lat2,lon2 = self.latitudeLastGridpoint,self.longitudeLastGridpoint
1463            if lon1 > 180.0: lon1 -= 360.0
1464            if lon2 > 180.0: lon2 -= 360.0
1465            lats = np.linspace(lat1,lat2,self.ny)
1466            lons = np.linspace(lon1,lon2,self.nx)
1467            lons,lats = np.meshgrid(lons,lats) # Make 2-d arrays.
1468            if unrotate:
1469                from grib2io.utils import rotated_grid
1470                lats,lons = rotated_grid.unrotate(lats,lons,self.anglePoleRotation,
1471                                                  self.latitudeSouthernPole,
1472                                                  self.longitudeSouthernPole)
1473        elif gdtn == 40: # Gaussian grid (only works for global!)
1474            from grib2io.utils.gauss_grid import gaussian_latitudes
1475            lon1, lat1 = self.longitudeFirstGridpoint, self.latitudeFirstGridpoint
1476            lon2, lat2 = self.longitudeLastGridpoint, self.latitudeLastGridpoint
1477            nlats = self.ny
1478            if not reggrid: # Reduced Gaussian grid.
1479                nlons = 2*nlats
1480                dlon = 360./nlons
1481            else:
1482                nlons = self.nx
1483                dlon = self.gridlengthXDirection
1484            lons = np.linspace(lon1,lon2,nlons)
1485            # Compute Gaussian lats (north to south)
1486            lats = gaussian_latitudes(nlats)
1487            if lat1 < lat2:  # reverse them if necessary
1488                lats = lats[::-1]
1489            lons,lats = np.meshgrid(lons,lats)
1490        elif gdtn in {10,20,30,31,110}:
1491            # Mercator, Lambert Conformal, Stereographic, Albers Equal Area,
1492            # Azimuthal Equidistant
1493            dx,dy = self.gridlengthXDirection, self.gridlengthYDirection
1494            lon1,lat1 = self.longitudeFirstGridpoint, self.latitudeFirstGridpoint
1495            pj = pyproj.Proj(self.projParameters)
1496            llcrnrx, llcrnry = pj(lon1,lat1)
1497            x = llcrnrx+dx*np.arange(self.nx)
1498            y = llcrnry+dy*np.arange(self.ny)
1499            x,y = np.meshgrid(x, y)
1500            lons,lats = pj(x, y, inverse=True)
1501        elif gdtn == 90:
1502            # Satellite Projection
1503            dx = self.gridlengthXDirection
1504            dy = self.gridlengthYDirection
1505            pj = pyproj.Proj(self.projParameters)
1506            x = dx*np.indices((self.ny,self.nx),'f')[1,:,:]
1507            x -= 0.5*x.max()
1508            y = dy*np.indices((self.ny,self.nx),'f')[0,:,:]
1509            y -= 0.5*y.max()
1510            lons,lats = pj(x,y,inverse=True)
1511            # Set lons,lats to 1.e30 where undefined
1512            abslons = np.fabs(lons)
1513            abslats = np.fabs(lats)
1514            lons = np.where(abslons < 1.e20, lons, 1.e30)
1515            lats = np.where(abslats < 1.e20, lats, 1.e30)
1516        elif gdtn == 32769:
1517            # Special NCEP Grid, Rotated Lat/Lon, Arakawa E-Grid (Non-Staggered)
1518            from grib2io.utils import arakawa_rotated_grid
1519            from grib2io.utils.rotated_grid import DEG2RAD
1520            di, dj = 0.0, 0.0
1521            do_180 = False
1522            idir = 1 if self.scanModeFlags[0] == 0 else -1
1523            jdir = -1 if self.scanModeFlags[1] == 0 else 1
1524            grid_oriented = 0 if self.resolutionAndComponentFlags[4] == 0 else 1
1525            do_rot = 1 if grid_oriented == 1 else 0
1526            la1 = self.latitudeFirstGridpoint
1527            lo1 = self.longitudeFirstGridpoint
1528            clon = self.longitudeCenterGridpoint
1529            clat = self.latitudeCenterGridpoint
1530            lasp = clat - 90.0
1531            losp = clon
1532            llat, llon = arakawa_rotated_grid.ll2rot(la1,lo1,lasp,losp)
1533            la2, lo2 = arakawa_rotated_grid.rot2ll(-llat,-llon,lasp,losp)
1534            rlat = -llat
1535            rlon = -llon
1536            if self.nx == 1:
1537                di = 0.0
1538            elif idir == 1:
1539                ti = rlon
1540                while ti < llon:
1541                    ti += 360.0
1542                di = (ti - llon)/float(self.nx-1)
1543            else:
1544                ti = llon
1545                while ti < rlon:
1546                    ti += 360.0
1547                di = (ti - rlon)/float(self.nx-1)
1548            if self.ny == 1:
1549               dj = 0.0
1550            else:
1551                dj = (rlat - llat)/float(self.ny-1)
1552                if dj < 0.0:
1553                    dj = -dj
1554            if idir == 1:
1555                if llon > rlon:
1556                    llon -= 360.0
1557                if llon < 0 and rlon > 0:
1558                    do_180 = True
1559            else:
1560                if rlon > llon:
1561                    rlon -= 360.0
1562                if rlon < 0 and llon > 0:
1563                    do_180 = True
1564            xlat1d = llat + (np.arange(self.ny)*jdir*dj)
1565            xlon1d = llon + (np.arange(self.nx)*idir*di)
1566            xlons, xlats = np.meshgrid(xlon1d,xlat1d)
1567            rot2ll_vectorized = np.vectorize(arakawa_rotated_grid.rot2ll)
1568            lats, lons = rot2ll_vectorized(xlats,xlons,lasp,losp)
1569            if do_180:
1570                lons = np.where(lons>180.0,lons-360.0,lons)
1571            vector_rotation_angles_vectorized = np.vectorize(arakawa_rotated_grid.vector_rotation_angles)
1572            rots = vector_rotation_angles_vectorized(lats, lons, clat, losp, xlats)
1573            del xlat1d, xlon1d, xlats, xlons
1574        else:
1575            raise ValueError('Unsupported grid')
1576
1577        _latlon_datastore[self._sha1_section3] = dict(latitude=lats,longitude=lons)
1578        try:
1579            _latlon_datastore[self._sha1_section3]['vector_rotation_angles'] = rots
1580        except(NameError):
1581            pass
1582
1583        return lats, lons
1584
1585
1586    def map_keys(self):
1587        """
1588        Unpack data grid replacing integer values with strings.
1589
1590        These types of fields are categorical or classifications where data
1591        values do not represent an observable or predictable physical quantity.
1592        An example of such a field would be [Dominant Precipitation Type -
1593        DPTYPE](https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_table4-201.shtml)
1594
1595        Returns
1596        -------
1597        map_keys
1598            numpy.ndarray of string values per element.
1599        """
1600        hold_auto_nans = _AUTO_NANS
1601        set_auto_nans(False)
1602
1603        if (np.all(self.section1[0:2] == [7, 14]) and self.shortName == "PWTHER") or (
1604            self._isNDFD and self.shortName in {"WX", "WWA"}
1605        ):
1606            keys = utils.decode_wx_strings(self.section2)
1607            if hasattr(self,'priMissingValue') and self.priMissingValue not in [None,0]:
1608                keys[int(self.priMissingValue)] = 'Missing'
1609            if hasattr(self,'secMissingValue') and self.secMissingValue not in [None,0]:
1610                keys[int(self.secMissingValue)] = 'Missing'
1611            u,inv = np.unique(self.data,return_inverse=True)
1612            fld = np.array([keys[x] for x in u])[inv].reshape(self.data.shape)
1613        else:
1614            # For data whose units are defined in a code table (i.e. classification or mask)
1615            tblname = re.findall(r'\d\.\d+',self.units,re.IGNORECASE)[0]
1616            fld = self.data.astype(np.int32).astype(str)
1617            tbl = tables.get_table(tblname,expand=True)
1618            for val in np.unique(fld):
1619                fld = np.where(fld==val,tbl[val],fld)
1620        set_auto_nans(hold_auto_nans)
1621        return fld
1622
1623
1624    def to_bytes(self, validate: bool=True):
1625        """
1626        Return packed GRIB2 message in bytes format.
1627
1628        This will be useful for exporting data in non-file formats. For example,
1629        can be used to output grib data directly to S3 using the boto3 client
1630        without the need to write a temporary file to upload first.
1631
1632        Parameters
1633        ----------
1634        validate: default=True
1635            If `True` (DEFAULT), validates first/last four bytes for proper
1636            formatting, else returns None. If `False`, message is output as is.
1637
1638        Returns
1639        -------
1640        to_bytes
1641            Returns GRIB2 formatted message as bytes.
1642        """
1643        if hasattr(self,'_msg'):
1644            if validate:
1645                if self.validate():
1646                    return self._msg
1647                else:
1648                    return None
1649            else:
1650                return self._msg
1651        else:
1652            return None
1653
1654
1655    def interpolate(self, method, grid_def_out, method_options=None, drtn=None,
1656                    num_threads=1):
1657        """
1658        Grib2Message Interpolator
1659
1660        Performs spatial interpolation via the [NCEPLIBS-ip
1661        library](https://github.com/NOAA-EMC/NCEPLIBS-ip). This interpolate
1662        method only supports scalar interpolation. If you need to perform
1663        vector interpolation, use the module-level `grib2io.interpolate`
1664        function.
1665
1666        Parameters
1667        ----------
1668        method
1669            Interpolate method to use. This can either be an integer or string
1670            using the following mapping:
1671
1672            | Interpolate Scheme | Integer Value |
1673            | :---:              | :---:         |
1674            | 'bilinear'         | 0             |
1675            | 'bicubic'          | 1             |
1676            | 'neighbor'         | 2             |
1677            | 'budget'           | 3             |
1678            | 'spectral'         | 4             |
1679            | 'neighbor-budget'  | 6             |
1680
1681        grid_def_out : grib2io.Grib2GridDef
1682            Grib2GridDef object of the output grid.
1683        method_options : list of ints, optional
1684            Interpolation options. See the NCEPLIBS-ip documentation for
1685            more information on how these are used.
1686        drtn
1687            Data Representation Template to be used for the returned
1688            interpolated GRIB2 message. When `None`, the data representation
1689            template of the source GRIB2 message is used. Once again, it is the
1690            user's responsibility to properly set the Data Representation
1691            Template attributes.
1692        num_threads : int, optional
1693            Number of OpenMP threads to use for interpolation. The default
1694            value is 1. If NCEPLIBS-ip and grib2io's iplib extension module
1695            was not built with OpenMP, then this keyword argument and value
1696            will have no impact.
1697
1698        Returns
1699        -------
1700        interpolate
1701            If interpolating to a grid, a new Grib2Message object is returned.
1702            The GRIB2 metadata of the new Grib2Message object is identical to
1703            the input except where required to be different because of the new
1704            grid specs and possibly a new data representation template.
1705
1706            If interpolating to station points, the interpolated data values are
1707            returned as a numpy.ndarray.
1708        """
1709        section0 = self.section0
1710        section0[-1] = 0
1711        gds = [0, grid_def_out.npoints, 0, 255, grid_def_out.gdtn]
1712        section3 = np.concatenate((gds,grid_def_out.gdt))
1713        drtn = self.drtn if drtn is None else drtn
1714
1715        msg = Grib2Message(section0,self.section1,self.section2,section3,
1716                           self.section4,None,self.bitMapFlag.value,drtn=drtn)
1717
1718        msg._msgnum = -1
1719        msg._deflist = self._deflist
1720        msg._coordlist = self._coordlist
1721        shape = (msg.ny,msg.nx)
1722        ndim = 2
1723        if msg.typeOfValues == 0:
1724            dtype = 'float32'
1725        elif msg.typeOfValues == 1:
1726            dtype = 'int32'
1727        msg._data = interpolate(self.data,method,Grib2GridDef.from_section3(self.section3),grid_def_out,
1728                                method_options=method_options,num_threads=num_threads).reshape(msg.ny,msg.nx)
1729        msg.section5[0] = grid_def_out.npoints
1730        return msg
1731
1732
1733    def subset(self, lats, lons):
1734        """
1735        Return a spatial subset.
1736
1737        Currently only supports regular grids of the following types:
1738
1739        | Grid Type                                                    | gdtn  |
1740        | :---:                                                        | :---: |
1741        | Latitude/Longitude, Equidistant Cylindrical, or Plate Carree | 0     |
1742        | Rotated Latitude/Longitude                                   | 1     |
1743        | Mercator                                                     | 10    |
1744        | Polar Stereographic                                          | 20    |
1745        | Lambert Conformal                                            | 30    |
1746        | Albers Equal-Area                                            | 31    |
1747        | Gaussian Latitude/Longitude                                  | 40    |
1748        | Equatorial Azimuthal Equidistant Projection                  | 110   |
1749
1750        Parameters
1751        ----------
1752        lats
1753            List or tuple of latitudes.  The minimum and maximum latitudes will
1754            be used to define the southern and northern boundaries.
1755
1756            The order of the latitudes is not important.  The function will
1757            determine which is the minimum and maximum.
1758
1759            The latitudes should be in decimal degrees with 0.0 at the equator,
1760            positive values in the northern hemisphere increasing to 90, and
1761            negative values in the southern hemisphere decreasing to -90.
1762        lons
1763            List or tuple of longitudes.  The minimum and maximum longitudes
1764            will be used to define the western and eastern boundaries.
1765
1766            The order of the longitudes is not important.  The function will
1767            determine which is the minimum and maximum.
1768
1769            GRIB2 longitudes should be in decimal degrees with 0.0 at the prime
1770            meridian, positive values increasing eastward to 360.  There are no
1771            negative GRIB2 longitudes.
1772
1773            The typical west longitudes that start at 0.0 at the prime meridian
1774            and decrease to -180 westward, are converted to GRIB2 longitudes by
1775            '360 - (absolute value of the west longitude)' where typical
1776            eastern longitudes are unchanged as GRIB2 longitudes.
1777
1778        Returns
1779        -------
1780        subset
1781            A spatial subset of a GRIB2 message.
1782        """
1783        if self.gdtn not in [0, 1, 10, 20, 30, 31, 40, 110]:
1784            raise ValueError(
1785                """
1786
1787Subset only works for
1788    Latitude/Longitude, Equidistant Cylindrical, or Plate Carree (gdtn=0)
1789    Rotated Latitude/Longitude (gdtn=1)
1790    Mercator (gdtn=10)
1791    Polar Stereographic (gdtn=20)
1792    Lambert Conformal (gdtn=30)
1793    Albers Equal-Area (gdtn=31)
1794    Gaussian Latitude/Longitude (gdtn=40)
1795    Equatorial Azimuthal Equidistant Projection (gdtn=110)
1796
1797"""
1798            )
1799
1800        if self.nx == 0 or self.ny == 0:
1801            raise ValueError(
1802                """
1803
1804Subset only works for regular grids.
1805
1806"""
1807            )
1808
1809        newmsg = Grib2Message(
1810            np.copy(self.section0),
1811            np.copy(self.section1),
1812            np.copy(self.section2),
1813            np.copy(self.section3),
1814            np.copy(self.section4),
1815            np.copy(self.section5),
1816        )
1817
1818        msglats, msglons = self.grid()
1819
1820        la1 = np.max(lats)
1821        lo1 = np.min(lons)
1822        la2 = np.min(lats)
1823        lo2 = np.max(lons)
1824
1825        # Find the indices of the first and last grid points to the nearest
1826        # lat/lon values in the grid.
1827        first_lat = np.abs(msglats - la1)
1828        first_lon = np.abs(msglons - lo1)
1829        max_idx = np.maximum(first_lat, first_lon)
1830        first_j, first_i = np.where(max_idx == np.min(max_idx))
1831
1832        last_lat = np.abs(msglats - la2)
1833        last_lon = np.abs(msglons - lo2)
1834        max_idx = np.maximum(last_lat, last_lon)
1835        last_j, last_i = np.where(max_idx == np.min(max_idx))
1836
1837        setattr(newmsg, "latitudeFirstGridpoint", msglats[first_j[0], first_i[0]])
1838        setattr(newmsg, "longitudeFirstGridpoint", msglons[first_j[0], first_i[0]])
1839        setattr(newmsg, "nx", np.abs(first_i[0] - last_i[0]))
1840        setattr(newmsg, "ny", np.abs(first_j[0] - last_j[0]))
1841
1842        # Set *LastGridpoint attributes even if only used for gdtn=[0, 1, 40].
1843        # This information is used to subset xarray datasets and even though
1844        # unnecessary for some supported grid types, it won't affect a grib2io
1845        # message to set them.
1846        setattr(newmsg, "latitudeLastGridpoint", msglats[last_j[0], last_i[0]])
1847        setattr(newmsg, "longitudeLastGridpoint", msglons[last_j[0], last_i[0]])
1848
1849        setattr(
1850            newmsg,
1851            "data",
1852            self.data[
1853                min(first_j[0], last_j[0]) : max(first_j[0], last_j[0]),
1854                min(first_i[0], last_i[0]) : max(first_i[0], last_i[0]),
1855            ].copy(),
1856        )
1857
1858        # Need to set the newmsg._sha1_section3 to a blank string so the grid
1859        # method ignores the cached lat/lon values.  This will force the grid
1860        # method to recompute the lat/lon values for the subsetted grid.
1861        newmsg._sha1_section3 = ""
1862        newmsg.grid()
1863
1864        return newmsg
1865
1866    def validate(self):
1867        """
1868        Validate a complete GRIB2 message.
1869
1870        The g2c library does its own internal validation when g2_gribend() is called, but
1871        we will check in grib2io also. The validation checks if the first 4 bytes in
1872        self._msg is 'GRIB' and '7777' as the last 4 bytes and that the message length in
1873        section 0 equals the length of the packed message.
1874
1875        Returns
1876        -------
1877        `True` if the packed GRIB2 message is complete and well-formed, `False` otherwise.
1878        """
1879        valid = False
1880        if hasattr(self,'_msg'):
1881            if self._msg[0:4]+self._msg[-4:] == b'GRIB7777':
1882                if self.section0[-1] == len(self._msg):
1883                    valid = True
1884        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
1074    @property
1075    def gdtn(self):
1076        """Return Grid Definition Template Number"""
1077        return self.section3[4]

Return Grid Definition Template Number

gdt
1080    @property
1081    def gdt(self):
1082        """Return Grid Definition Template."""
1083        return self.gridDefinitionTemplate

Return Grid Definition Template.

pdtn
1086    @property
1087    def pdtn(self):
1088        """Return Product Definition Template Number."""
1089        return self.section4[1]

Return Product Definition Template Number.

pdt
1092    @property
1093    def pdt(self):
1094        """Return Product Definition Template."""
1095        return self.productDefinitionTemplate

Return Product Definition Template.

drtn
1098    @property
1099    def drtn(self):
1100        """Return Data Representation Template Number."""
1101        return self.section5[1]

Return Data Representation Template Number.

drt
1104    @property
1105    def drt(self):
1106        """Return Data Representation Template."""
1107        return self.dataRepresentationTemplate

Return Data Representation Template.

pdy
1110    @property
1111    def pdy(self):
1112        """Return the PDY ('YYYYMMDD')."""
1113        return ''.join([str(i) for i in self.section1[5:8]])

Return the PDY ('YYYYMMDD').

griddef
1116    @property
1117    def griddef(self):
1118        """Return a Grib2GridDef instance for a GRIB2 message."""
1119        return Grib2GridDef.from_section3(self.section3)

Return a Grib2GridDef instance for a GRIB2 message.

lats
1122    @property
1123    def lats(self):
1124        """Return grid latitudes."""
1125        return self.latlons()[0]

Return grid latitudes.

lons
1128    @property
1129    def lons(self):
1130        """Return grid longitudes."""
1131        return self.latlons()[1]

Return grid longitudes.

min
1133    @property
1134    def min(self):
1135        """Return minimum value of data."""
1136        return np.nanmin(self.data)

Return minimum value of data.

max
1139    @property
1140    def max(self):
1141        """Return maximum value of data."""
1142        return np.nanmax(self.data)

Return maximum value of data.

mean
1145    @property
1146    def mean(self):
1147        """Return mean value of data."""
1148        return np.nanmean(self.data)

Return mean value of data.

median
1151    @property
1152    def median(self):
1153        """Return median value of data."""
1154        return np.nanmedian(self.data)

Return median value of data.

shape
1157    @property
1158    def shape(self):
1159        """Return shape of data."""
1160        return self.griddef.shape

Return shape of data.

def attrs_by_section(self, sect: int, values: bool = False):
1201    def attrs_by_section(self, sect: int, values: bool=False):
1202        """
1203        Provide a tuple of attribute names for the given GRIB2 section.
1204
1205        Parameters
1206        ----------
1207        sect
1208            The GRIB2 section number.
1209        values
1210            Optional (default is `False`) argument to return attributes values.
1211
1212        Returns
1213        -------
1214        attrs_by_section
1215            A list of attribute names or dict of name:value pairs if `values =
1216            True`.
1217        """
1218        if sect in {0,1,6}:
1219            attrs = templates._section_attrs[sect]
1220        elif sect in {3,4,5}:
1221            def _find_class_index(n):
1222                _key = {3:'Grid', 4:'Product', 5:'Data'}
1223                for i,c in enumerate(self.__class__.__mro__):
1224                    if _key[n] in c.__name__:
1225                        return i
1226                else:
1227                    return []
1228            if sys.version_info.minor <= 8:
1229                attrs = templates._section_attrs[sect]+\
1230                        [a for a in dir(self.__class__.__mro__[_find_class_index(sect)]) if not a.startswith('_')]
1231            else:
1232                attrs = templates._section_attrs[sect]+\
1233                        self.__class__.__mro__[_find_class_index(sect)]._attrs()
1234        else:
1235            attrs = []
1236        if values:
1237            return {k:getattr(self,k) for k in attrs}
1238        else:
1239            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):
1242    def copy(self, deep: bool = True):
1243        """Returns a copy of this Grib2Message.
1244
1245        When `deep=True`, a copy is made of each of the GRIB2 section arrays and
1246        the data are unpacked from the source object and copied into the new
1247        object. Otherwise, a shallow copy of each array is performed and no data
1248        are copied.
1249
1250        Parameters
1251        ----------
1252        deep : bool, default: True
1253            Whether each GRIB2 section array and data are copied onto
1254            the new object. Default is True.
1255
1256        Returns
1257        -------
1258        object : Grib2Message
1259            New Grib2Message object.
1260
1261            .. versionadded:: 2.6.0
1262        """
1263        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):
1266    def pack(self):
1267        """
1268        Pack GRIB2 section data into a binary message.
1269
1270        It is the user's responsibility to populate the GRIB2 section
1271        information with appropriate metadata.
1272        """
1273        # Create beginning of packed binary message with section 0 and 1 data.
1274        self._sections = []
1275        self._msg,self._pos = g2clib.grib2_create(self.indicatorSection[2:4],self.identificationSection)
1276        self._sections += [0,1]
1277
1278        # Add section 2 if present.
1279        if isinstance(self.section2,bytes) and len(self.section2) > 0:
1280            self._msg,self._pos = g2clib.grib2_addlocal(self._msg,self.section2)
1281            self._sections.append(2)
1282
1283        # Add section 3.
1284        self.section3[1] = self.nx * self.ny
1285        self._msg,self._pos = g2clib.grib2_addgrid(self._msg,self.gridDefinitionSection,
1286                                                   self.gridDefinitionTemplate,
1287                                                   self._deflist)
1288        self._sections.append(3)
1289
1290        # Prepare data.
1291        if self._data is None:
1292            if self._ondiskarray is None:
1293                raise ValueError("Grib2Message object has no data, thus it cannot be packed.")
1294        field = np.copy(self.data)
1295        if self.scanModeFlags is not None:
1296            if self.scanModeFlags[3]:
1297                fieldsave = field.astype('f') # Casting makes a copy
1298                field[1::2,:] = fieldsave[1::2,::-1]
1299        fld = field.astype('f')
1300
1301        # Prepare bitmap, if necessary
1302        bitmapflag = self.bitMapFlag.value
1303        if bitmapflag == 0:
1304            if self.bitmap is not None:
1305                bmap = np.ravel(self.bitmap).astype(DEFAULT_NUMPY_INT)
1306            else:
1307                bmap = np.ravel(np.where(np.isnan(fld),0,1)).astype(DEFAULT_NUMPY_INT)
1308        else:
1309            bmap = None
1310
1311        # Prepare data for packing if nans are present
1312        fld = np.ravel(fld)
1313        if bitmapflag in {0,254}:
1314            fld = np.where(np.isnan(fld),0,fld)
1315        else:
1316            if np.isnan(fld).any():
1317                if hasattr(self,'priMissingValue'):
1318                    fld = np.where(np.isnan(fld),self.priMissingValue,fld)
1319            if hasattr(self,'_missvalmap'):
1320                if hasattr(self,'priMissingValue'):
1321                    fld = np.where(self._missvalmap==1,self.priMissingValue,fld)
1322                if hasattr(self,'secMissingValue'):
1323                    fld = np.where(self._missvalmap==2,self.secMissingValue,fld)
1324
1325        # Add sections 4, 5, 6, and 7.
1326        self._msg,self._pos = g2clib.grib2_addfield(self._msg,self.pdtn,
1327                                                    self.productDefinitionTemplate,
1328                                                    self._coordlist,
1329                                                    self.drtn,
1330                                                    self.dataRepresentationTemplate,
1331                                                    fld,
1332                                                    bitmapflag,
1333                                                    bmap)
1334        self._sections.append(4)
1335        self._sections.append(5)
1336        self._sections.append(6)
1337        self._sections.append(7)
1338
1339        # Finalize GRIB2 message with section 8.
1340        self._msg, self._pos = g2clib.grib2_end(self._msg)
1341        self._sections.append(8)
1342        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>
1345    @property
1346    def data(self) -> np.array:
1347        """Access the unpacked data values."""
1348        if self._data is None:
1349            if self._auto_nans != _AUTO_NANS:
1350                self._data = self._ondiskarray
1351            self._data = np.asarray(self._ondiskarray)
1352        return self._data

Access the unpacked data values.

def flush_data(self):
1393    def flush_data(self):
1394        """
1395        Flush the unpacked data values from the Grib2Message object.
1396
1397        Notes
1398        -----
1399        If the Grib2Message object was constructed from "scratch" (i.e.
1400        not read from file), this method will remove the data array from
1401        the object and it cannot be recovered.
1402        """
1403        self._data = None
1404        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):
1415    def latlons(self, *args, **kwrgs):
1416        """Alias for `grib2io.Grib2Message.grid` method."""
1417        return self.grid(*args, **kwrgs)

Alias for grib2io.Grib2Message.grid method.

def grid(self, unrotate: bool = True):
1420    def grid(self, unrotate: bool=True):
1421        """
1422        Return lats,lons (in degrees) of grid.
1423
1424        Currently can handle reg. lat/lon,cglobal Gaussian, mercator,
1425        stereographic, lambert conformal, albers equal-area, space-view and
1426        azimuthal equidistant grids.
1427
1428        Parameters
1429        ----------
1430        unrotate
1431            If `True` [DEFAULT], and grid is rotated lat/lon, then unrotate the
1432            grid, otherwise `False`, do not.
1433
1434        Returns
1435        -------
1436        lats, lons : numpy.ndarray
1437            Returns two numpy.ndarrays with dtype=numpy.float32 of grid
1438            latitudes and longitudes in units of degrees.
1439        """
1440        if self._sha1_section3 in _latlon_datastore.keys():
1441            return (_latlon_datastore[self._sha1_section3]['latitude'],
1442                    _latlon_datastore[self._sha1_section3]['longitude'])
1443        gdtn = self.gridDefinitionTemplateNumber.value
1444        gdtmpl = self.gridDefinitionTemplate
1445        reggrid = self.gridDefinitionSection[2] == 0 # This means regular 2-d grid
1446        if gdtn == 0:
1447            # Regular lat/lon grid
1448            lon1, lat1 = self.longitudeFirstGridpoint, self.latitudeFirstGridpoint
1449            lon2, lat2 = self.longitudeLastGridpoint, self.latitudeLastGridpoint
1450            dlon = self.gridlengthXDirection
1451            dlat = self.gridlengthYDirection
1452            if lon2 < lon1 and dlon < 0: lon1 = -lon1
1453            lats = np.linspace(lat1,lat2,self.ny)
1454            if reggrid:
1455                lons = np.linspace(lon1,lon2,self.nx)
1456            else:
1457                lons = np.linspace(lon1,lon2,self.ny*2)
1458            lons,lats = np.meshgrid(lons,lats) # Make 2-d arrays.
1459        elif gdtn == 1: # Rotated Lat/Lon grid
1460            pj = pyproj.Proj(self.projParameters)
1461            lat1,lon1 = self.latitudeFirstGridpoint,self.longitudeFirstGridpoint
1462            lat2,lon2 = self.latitudeLastGridpoint,self.longitudeLastGridpoint
1463            if lon1 > 180.0: lon1 -= 360.0
1464            if lon2 > 180.0: lon2 -= 360.0
1465            lats = np.linspace(lat1,lat2,self.ny)
1466            lons = np.linspace(lon1,lon2,self.nx)
1467            lons,lats = np.meshgrid(lons,lats) # Make 2-d arrays.
1468            if unrotate:
1469                from grib2io.utils import rotated_grid
1470                lats,lons = rotated_grid.unrotate(lats,lons,self.anglePoleRotation,
1471                                                  self.latitudeSouthernPole,
1472                                                  self.longitudeSouthernPole)
1473        elif gdtn == 40: # Gaussian grid (only works for global!)
1474            from grib2io.utils.gauss_grid import gaussian_latitudes
1475            lon1, lat1 = self.longitudeFirstGridpoint, self.latitudeFirstGridpoint
1476            lon2, lat2 = self.longitudeLastGridpoint, self.latitudeLastGridpoint
1477            nlats = self.ny
1478            if not reggrid: # Reduced Gaussian grid.
1479                nlons = 2*nlats
1480                dlon = 360./nlons
1481            else:
1482                nlons = self.nx
1483                dlon = self.gridlengthXDirection
1484            lons = np.linspace(lon1,lon2,nlons)
1485            # Compute Gaussian lats (north to south)
1486            lats = gaussian_latitudes(nlats)
1487            if lat1 < lat2:  # reverse them if necessary
1488                lats = lats[::-1]
1489            lons,lats = np.meshgrid(lons,lats)
1490        elif gdtn in {10,20,30,31,110}:
1491            # Mercator, Lambert Conformal, Stereographic, Albers Equal Area,
1492            # Azimuthal Equidistant
1493            dx,dy = self.gridlengthXDirection, self.gridlengthYDirection
1494            lon1,lat1 = self.longitudeFirstGridpoint, self.latitudeFirstGridpoint
1495            pj = pyproj.Proj(self.projParameters)
1496            llcrnrx, llcrnry = pj(lon1,lat1)
1497            x = llcrnrx+dx*np.arange(self.nx)
1498            y = llcrnry+dy*np.arange(self.ny)
1499            x,y = np.meshgrid(x, y)
1500            lons,lats = pj(x, y, inverse=True)
1501        elif gdtn == 90:
1502            # Satellite Projection
1503            dx = self.gridlengthXDirection
1504            dy = self.gridlengthYDirection
1505            pj = pyproj.Proj(self.projParameters)
1506            x = dx*np.indices((self.ny,self.nx),'f')[1,:,:]
1507            x -= 0.5*x.max()
1508            y = dy*np.indices((self.ny,self.nx),'f')[0,:,:]
1509            y -= 0.5*y.max()
1510            lons,lats = pj(x,y,inverse=True)
1511            # Set lons,lats to 1.e30 where undefined
1512            abslons = np.fabs(lons)
1513            abslats = np.fabs(lats)
1514            lons = np.where(abslons < 1.e20, lons, 1.e30)
1515            lats = np.where(abslats < 1.e20, lats, 1.e30)
1516        elif gdtn == 32769:
1517            # Special NCEP Grid, Rotated Lat/Lon, Arakawa E-Grid (Non-Staggered)
1518            from grib2io.utils import arakawa_rotated_grid
1519            from grib2io.utils.rotated_grid import DEG2RAD
1520            di, dj = 0.0, 0.0
1521            do_180 = False
1522            idir = 1 if self.scanModeFlags[0] == 0 else -1
1523            jdir = -1 if self.scanModeFlags[1] == 0 else 1
1524            grid_oriented = 0 if self.resolutionAndComponentFlags[4] == 0 else 1
1525            do_rot = 1 if grid_oriented == 1 else 0
1526            la1 = self.latitudeFirstGridpoint
1527            lo1 = self.longitudeFirstGridpoint
1528            clon = self.longitudeCenterGridpoint
1529            clat = self.latitudeCenterGridpoint
1530            lasp = clat - 90.0
1531            losp = clon
1532            llat, llon = arakawa_rotated_grid.ll2rot(la1,lo1,lasp,losp)
1533            la2, lo2 = arakawa_rotated_grid.rot2ll(-llat,-llon,lasp,losp)
1534            rlat = -llat
1535            rlon = -llon
1536            if self.nx == 1:
1537                di = 0.0
1538            elif idir == 1:
1539                ti = rlon
1540                while ti < llon:
1541                    ti += 360.0
1542                di = (ti - llon)/float(self.nx-1)
1543            else:
1544                ti = llon
1545                while ti < rlon:
1546                    ti += 360.0
1547                di = (ti - rlon)/float(self.nx-1)
1548            if self.ny == 1:
1549               dj = 0.0
1550            else:
1551                dj = (rlat - llat)/float(self.ny-1)
1552                if dj < 0.0:
1553                    dj = -dj
1554            if idir == 1:
1555                if llon > rlon:
1556                    llon -= 360.0
1557                if llon < 0 and rlon > 0:
1558                    do_180 = True
1559            else:
1560                if rlon > llon:
1561                    rlon -= 360.0
1562                if rlon < 0 and llon > 0:
1563                    do_180 = True
1564            xlat1d = llat + (np.arange(self.ny)*jdir*dj)
1565            xlon1d = llon + (np.arange(self.nx)*idir*di)
1566            xlons, xlats = np.meshgrid(xlon1d,xlat1d)
1567            rot2ll_vectorized = np.vectorize(arakawa_rotated_grid.rot2ll)
1568            lats, lons = rot2ll_vectorized(xlats,xlons,lasp,losp)
1569            if do_180:
1570                lons = np.where(lons>180.0,lons-360.0,lons)
1571            vector_rotation_angles_vectorized = np.vectorize(arakawa_rotated_grid.vector_rotation_angles)
1572            rots = vector_rotation_angles_vectorized(lats, lons, clat, losp, xlats)
1573            del xlat1d, xlon1d, xlats, xlons
1574        else:
1575            raise ValueError('Unsupported grid')
1576
1577        _latlon_datastore[self._sha1_section3] = dict(latitude=lats,longitude=lons)
1578        try:
1579            _latlon_datastore[self._sha1_section3]['vector_rotation_angles'] = rots
1580        except(NameError):
1581            pass
1582
1583        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):
1586    def map_keys(self):
1587        """
1588        Unpack data grid replacing integer values with strings.
1589
1590        These types of fields are categorical or classifications where data
1591        values do not represent an observable or predictable physical quantity.
1592        An example of such a field would be [Dominant Precipitation Type -
1593        DPTYPE](https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_table4-201.shtml)
1594
1595        Returns
1596        -------
1597        map_keys
1598            numpy.ndarray of string values per element.
1599        """
1600        hold_auto_nans = _AUTO_NANS
1601        set_auto_nans(False)
1602
1603        if (np.all(self.section1[0:2] == [7, 14]) and self.shortName == "PWTHER") or (
1604            self._isNDFD and self.shortName in {"WX", "WWA"}
1605        ):
1606            keys = utils.decode_wx_strings(self.section2)
1607            if hasattr(self,'priMissingValue') and self.priMissingValue not in [None,0]:
1608                keys[int(self.priMissingValue)] = 'Missing'
1609            if hasattr(self,'secMissingValue') and self.secMissingValue not in [None,0]:
1610                keys[int(self.secMissingValue)] = 'Missing'
1611            u,inv = np.unique(self.data,return_inverse=True)
1612            fld = np.array([keys[x] for x in u])[inv].reshape(self.data.shape)
1613        else:
1614            # For data whose units are defined in a code table (i.e. classification or mask)
1615            tblname = re.findall(r'\d\.\d+',self.units,re.IGNORECASE)[0]
1616            fld = self.data.astype(np.int32).astype(str)
1617            tbl = tables.get_table(tblname,expand=True)
1618            for val in np.unique(fld):
1619                fld = np.where(fld==val,tbl[val],fld)
1620        set_auto_nans(hold_auto_nans)
1621        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):
1624    def to_bytes(self, validate: bool=True):
1625        """
1626        Return packed GRIB2 message in bytes format.
1627
1628        This will be useful for exporting data in non-file formats. For example,
1629        can be used to output grib data directly to S3 using the boto3 client
1630        without the need to write a temporary file to upload first.
1631
1632        Parameters
1633        ----------
1634        validate: default=True
1635            If `True` (DEFAULT), validates first/last four bytes for proper
1636            formatting, else returns None. If `False`, message is output as is.
1637
1638        Returns
1639        -------
1640        to_bytes
1641            Returns GRIB2 formatted message as bytes.
1642        """
1643        if hasattr(self,'_msg'):
1644            if validate:
1645                if self.validate():
1646                    return self._msg
1647                else:
1648                    return None
1649            else:
1650                return self._msg
1651        else:
1652            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):
1655    def interpolate(self, method, grid_def_out, method_options=None, drtn=None,
1656                    num_threads=1):
1657        """
1658        Grib2Message Interpolator
1659
1660        Performs spatial interpolation via the [NCEPLIBS-ip
1661        library](https://github.com/NOAA-EMC/NCEPLIBS-ip). This interpolate
1662        method only supports scalar interpolation. If you need to perform
1663        vector interpolation, use the module-level `grib2io.interpolate`
1664        function.
1665
1666        Parameters
1667        ----------
1668        method
1669            Interpolate method to use. This can either be an integer or string
1670            using the following mapping:
1671
1672            | Interpolate Scheme | Integer Value |
1673            | :---:              | :---:         |
1674            | 'bilinear'         | 0             |
1675            | 'bicubic'          | 1             |
1676            | 'neighbor'         | 2             |
1677            | 'budget'           | 3             |
1678            | 'spectral'         | 4             |
1679            | 'neighbor-budget'  | 6             |
1680
1681        grid_def_out : grib2io.Grib2GridDef
1682            Grib2GridDef object of the output grid.
1683        method_options : list of ints, optional
1684            Interpolation options. See the NCEPLIBS-ip documentation for
1685            more information on how these are used.
1686        drtn
1687            Data Representation Template to be used for the returned
1688            interpolated GRIB2 message. When `None`, the data representation
1689            template of the source GRIB2 message is used. Once again, it is the
1690            user's responsibility to properly set the Data Representation
1691            Template attributes.
1692        num_threads : int, optional
1693            Number of OpenMP threads to use for interpolation. The default
1694            value is 1. If NCEPLIBS-ip and grib2io's iplib extension module
1695            was not built with OpenMP, then this keyword argument and value
1696            will have no impact.
1697
1698        Returns
1699        -------
1700        interpolate
1701            If interpolating to a grid, a new Grib2Message object is returned.
1702            The GRIB2 metadata of the new Grib2Message object is identical to
1703            the input except where required to be different because of the new
1704            grid specs and possibly a new data representation template.
1705
1706            If interpolating to station points, the interpolated data values are
1707            returned as a numpy.ndarray.
1708        """
1709        section0 = self.section0
1710        section0[-1] = 0
1711        gds = [0, grid_def_out.npoints, 0, 255, grid_def_out.gdtn]
1712        section3 = np.concatenate((gds,grid_def_out.gdt))
1713        drtn = self.drtn if drtn is None else drtn
1714
1715        msg = Grib2Message(section0,self.section1,self.section2,section3,
1716                           self.section4,None,self.bitMapFlag.value,drtn=drtn)
1717
1718        msg._msgnum = -1
1719        msg._deflist = self._deflist
1720        msg._coordlist = self._coordlist
1721        shape = (msg.ny,msg.nx)
1722        ndim = 2
1723        if msg.typeOfValues == 0:
1724            dtype = 'float32'
1725        elif msg.typeOfValues == 1:
1726            dtype = 'int32'
1727        msg._data = interpolate(self.data,method,Grib2GridDef.from_section3(self.section3),grid_def_out,
1728                                method_options=method_options,num_threads=num_threads).reshape(msg.ny,msg.nx)
1729        msg.section5[0] = grid_def_out.npoints
1730        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):
1733    def subset(self, lats, lons):
1734        """
1735        Return a spatial subset.
1736
1737        Currently only supports regular grids of the following types:
1738
1739        | Grid Type                                                    | gdtn  |
1740        | :---:                                                        | :---: |
1741        | Latitude/Longitude, Equidistant Cylindrical, or Plate Carree | 0     |
1742        | Rotated Latitude/Longitude                                   | 1     |
1743        | Mercator                                                     | 10    |
1744        | Polar Stereographic                                          | 20    |
1745        | Lambert Conformal                                            | 30    |
1746        | Albers Equal-Area                                            | 31    |
1747        | Gaussian Latitude/Longitude                                  | 40    |
1748        | Equatorial Azimuthal Equidistant Projection                  | 110   |
1749
1750        Parameters
1751        ----------
1752        lats
1753            List or tuple of latitudes.  The minimum and maximum latitudes will
1754            be used to define the southern and northern boundaries.
1755
1756            The order of the latitudes is not important.  The function will
1757            determine which is the minimum and maximum.
1758
1759            The latitudes should be in decimal degrees with 0.0 at the equator,
1760            positive values in the northern hemisphere increasing to 90, and
1761            negative values in the southern hemisphere decreasing to -90.
1762        lons
1763            List or tuple of longitudes.  The minimum and maximum longitudes
1764            will be used to define the western and eastern boundaries.
1765
1766            The order of the longitudes is not important.  The function will
1767            determine which is the minimum and maximum.
1768
1769            GRIB2 longitudes should be in decimal degrees with 0.0 at the prime
1770            meridian, positive values increasing eastward to 360.  There are no
1771            negative GRIB2 longitudes.
1772
1773            The typical west longitudes that start at 0.0 at the prime meridian
1774            and decrease to -180 westward, are converted to GRIB2 longitudes by
1775            '360 - (absolute value of the west longitude)' where typical
1776            eastern longitudes are unchanged as GRIB2 longitudes.
1777
1778        Returns
1779        -------
1780        subset
1781            A spatial subset of a GRIB2 message.
1782        """
1783        if self.gdtn not in [0, 1, 10, 20, 30, 31, 40, 110]:
1784            raise ValueError(
1785                """
1786
1787Subset only works for
1788    Latitude/Longitude, Equidistant Cylindrical, or Plate Carree (gdtn=0)
1789    Rotated Latitude/Longitude (gdtn=1)
1790    Mercator (gdtn=10)
1791    Polar Stereographic (gdtn=20)
1792    Lambert Conformal (gdtn=30)
1793    Albers Equal-Area (gdtn=31)
1794    Gaussian Latitude/Longitude (gdtn=40)
1795    Equatorial Azimuthal Equidistant Projection (gdtn=110)
1796
1797"""
1798            )
1799
1800        if self.nx == 0 or self.ny == 0:
1801            raise ValueError(
1802                """
1803
1804Subset only works for regular grids.
1805
1806"""
1807            )
1808
1809        newmsg = Grib2Message(
1810            np.copy(self.section0),
1811            np.copy(self.section1),
1812            np.copy(self.section2),
1813            np.copy(self.section3),
1814            np.copy(self.section4),
1815            np.copy(self.section5),
1816        )
1817
1818        msglats, msglons = self.grid()
1819
1820        la1 = np.max(lats)
1821        lo1 = np.min(lons)
1822        la2 = np.min(lats)
1823        lo2 = np.max(lons)
1824
1825        # Find the indices of the first and last grid points to the nearest
1826        # lat/lon values in the grid.
1827        first_lat = np.abs(msglats - la1)
1828        first_lon = np.abs(msglons - lo1)
1829        max_idx = np.maximum(first_lat, first_lon)
1830        first_j, first_i = np.where(max_idx == np.min(max_idx))
1831
1832        last_lat = np.abs(msglats - la2)
1833        last_lon = np.abs(msglons - lo2)
1834        max_idx = np.maximum(last_lat, last_lon)
1835        last_j, last_i = np.where(max_idx == np.min(max_idx))
1836
1837        setattr(newmsg, "latitudeFirstGridpoint", msglats[first_j[0], first_i[0]])
1838        setattr(newmsg, "longitudeFirstGridpoint", msglons[first_j[0], first_i[0]])
1839        setattr(newmsg, "nx", np.abs(first_i[0] - last_i[0]))
1840        setattr(newmsg, "ny", np.abs(first_j[0] - last_j[0]))
1841
1842        # Set *LastGridpoint attributes even if only used for gdtn=[0, 1, 40].
1843        # This information is used to subset xarray datasets and even though
1844        # unnecessary for some supported grid types, it won't affect a grib2io
1845        # message to set them.
1846        setattr(newmsg, "latitudeLastGridpoint", msglats[last_j[0], last_i[0]])
1847        setattr(newmsg, "longitudeLastGridpoint", msglons[last_j[0], last_i[0]])
1848
1849        setattr(
1850            newmsg,
1851            "data",
1852            self.data[
1853                min(first_j[0], last_j[0]) : max(first_j[0], last_j[0]),
1854                min(first_i[0], last_i[0]) : max(first_i[0], last_i[0]),
1855            ].copy(),
1856        )
1857
1858        # Need to set the newmsg._sha1_section3 to a blank string so the grid
1859        # method ignores the cached lat/lon values.  This will force the grid
1860        # method to recompute the lat/lon values for the subsetted grid.
1861        newmsg._sha1_section3 = ""
1862        newmsg.grid()
1863
1864        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):
1866    def validate(self):
1867        """
1868        Validate a complete GRIB2 message.
1869
1870        The g2c library does its own internal validation when g2_gribend() is called, but
1871        we will check in grib2io also. The validation checks if the first 4 bytes in
1872        self._msg is 'GRIB' and '7777' as the last 4 bytes and that the message length in
1873        section 0 equals the length of the packed message.
1874
1875        Returns
1876        -------
1877        `True` if the packed GRIB2 message is complete and well-formed, `False` otherwise.
1878        """
1879        valid = False
1880        if hasattr(self,'_msg'):
1881            if self._msg[0:4]+self._msg[-4:] == b'GRIB7777':
1882                if self.section0[-1] == len(self._msg):
1883                    valid = True
1884        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:
2381@dataclass
2382class Grib2GridDef:
2383    """
2384    Class for Grid Definition Template Number and Template as attributes.
2385
2386    This allows for cleaner looking code when passing these metadata around.
2387    For example, the `grib2io._Grib2Message.interpolate` method and
2388    `grib2io.interpolate` function accepts these objects.
2389    """
2390    gdtn: int
2391    gdt: NDArray
2392
2393    @classmethod
2394    def from_section3(cls, section3):
2395        return cls(section3[4],section3[5:])
2396
2397    @property
2398    def nx(self):
2399        """Number of grid points in x-direction."""
2400        return int(self.gdt[7])
2401
2402    @property
2403    def ny(self):
2404        """Number of grid points in y-direction."""
2405        return int(self.gdt[8])
2406
2407    @property
2408    def npoints(self):
2409        """Total number of grid points."""
2410        return int(self.gdt[7] * self.gdt[8])
2411
2412    @property
2413    def shape(self):
2414        """Shape of the grid."""
2415        return (int(self.ny), int(self.nx))
2416
2417    def to_section3(self):
2418        """Return a full GRIB2 section3 array."""
2419        return np.array(
2420            [0, self.npoints, 0, 0, self.gdtn] + list(self.gdt)
2421        ).astype(np.int64)

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):
2393    @classmethod
2394    def from_section3(cls, section3):
2395        return cls(section3[4],section3[5:])
nx
2397    @property
2398    def nx(self):
2399        """Number of grid points in x-direction."""
2400        return int(self.gdt[7])

Number of grid points in x-direction.

ny
2402    @property
2403    def ny(self):
2404        """Number of grid points in y-direction."""
2405        return int(self.gdt[8])

Number of grid points in y-direction.

npoints
2407    @property
2408    def npoints(self):
2409        """Total number of grid points."""
2410        return int(self.gdt[7] * self.gdt[8])

Total number of grid points.

shape
2412    @property
2413    def shape(self):
2414        """Shape of the grid."""
2415        return (int(self.ny), int(self.nx))

Shape of the grid.

def to_section3(self):
2417    def to_section3(self):
2418        """Return a full GRIB2 section3 array."""
2419        return np.array(
2420            [0, self.npoints, 0, 0, self.gdtn] + list(self.gdt)
2421        ).astype(np.int64)

Return a full GRIB2 section3 array.

def msgs_from_index(index: dict, filehandle=None):
774def msgs_from_index(index: dict, filehandle=None):
775    """
776    Construct a list of Grib2Message objects from an index dictionary.
777
778    This function reconstructs a sequence of `Grib2Message` instances using
779    metadata sections stored in an index dictionary. If an open file handle is
780    provided, each message is linked to its on-disk binary data through a
781    `Grib2MessageOnDiskArray`, allowing deferred reading of the actual data
782    values from the GRIB2 file.
783
784    Parameters
785    ----------
786    index : dict
787        Dictionary containing parsed GRIB2 index information, including
788        section data arrays such as ``section0`` through ``section5``,
789        ``sectionOffset``, ``offset``, and ``bmapflag``.
790    filehandle : file-like object, optional
791        An open binary file handle to the GRIB2 file corresponding to the index.
792        If provided, the returned messages can access on-disk data arrays via
793        memory offsets. If not provided, only metadata will be available.
794
795    Returns
796    -------
797    list of Grib2Message
798        List of reconstructed `Grib2Message` objects built from the provided
799        index. Each message contains metadata, and if `filehandle` is given,
800        also references to on-disk data through a `Grib2MessageOnDiskArray`.
801
802    Notes
803    -----
804    - Each message is constructed by zipping the corresponding section entries
805      (sections 0–5 and bitmap flags).
806    - When a file handle is supplied, each message’s `_ondiskarray` attribute is
807      initialized to allow direct access to GRIB2 data values without loading
808      them fully into memory.
809    - The `_msgnum` attribute of each message is assigned sequentially to
810      preserve message order.
811    """
812    zipped = zip(
813        index["section0"],
814        index["section1"],
815        index["section2"],
816        index["section3"],
817        index["section4"],
818        index["section5"],
819        index["bmapflag"]
820    )
821    msgs = [Grib2Message(*sections) for sections in zipped]
822
823    if filehandle is not None:
824        for n, (msg, offset, secpos) in enumerate(zip(msgs, index["offset"], index["sectionOffset"])):
825            msg._ondiskarray = Grib2MessageOnDiskArray(
826                shape=(msg.ny, msg.nx),
827                ndim=2,
828                dtype=TYPE_OF_VALUES_DTYPE[msg.typeOfValues],
829                filehandle=filehandle,
830                msg=msg,
831                offset=offset,
832                bitmap_offset=secpos[6],
833                data_offset=secpos[7]
834            )
835            msg._msgnum = n
836    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.