Versions Compared

Key

  • This line was added.
  • This line was removed.
  • Formatting was changed.

...

Tabs Container
directionhorizontal


Tabs Page
titleFortran 90


Code Block
titlebufr_read_tempf.f90
! (C) Copyright 2005- ECMWF.
!
! This software is licensed under the terms of the Apache Licence Version 2.0
! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
!
! In applying this licence, ECMWF does not waive the privileges and immunities granted to it by
! virtue of its status as an intergovernmental organisation nor does it submit to any jurisdiction.
!
!
! Description: read and print radiosonde data from TEMP BUFR messages.
! If available this version also lists the position information from the WMO list
! (now OSCAR/Surface) appended to the reports by ECMWF
!
! Author: Bruce Ingleby
!
! Please note that TEMP reports can be encoded in various ways in BUFR. Therefore the code
! below might not work directly for other types of TEMP messages than the one used in the
! example. It is advised to use bufr_dump first to understand the structure of these messages.
!
program bufr_read_tempf
  use eccodes
  implicit none
  integer            :: ifile
  integer            :: iret
  integer            :: ibufr
  integer            :: i, count = 0
  integer            :: iflag
  integer            :: status_id, status_ht, status_time = 0, status_p
  integer            :: status_rsno, status_rssoft, status_balloonwt, statid_missing
  integer(kind=4)    :: sizews
  integer(kind=4)    :: blockNumber, stationNumber
  integer(kind=4)    :: ymd, hms
  logical            :: llstdonly = .True.  ! Set True to list standard levels only
  logical            :: llskip
  real(kind=8)       :: year, month, day, hour, minute, second
  real(kind=8)       :: htg, htp, htec = 0, sondeType, balloonwt
  ! real(kind=8), dimension(:), allocatable :: descriptors
  real(kind=8), dimension(:), allocatable :: lat, lon
  real(kind=8), dimension(:), allocatable :: timeVal, dlatVal, dlonVal, vssVal
  real(kind=8), dimension(:), allocatable :: presVal, zVal, tVal, tdVal, wdirVal, wspVal
  character(len=128) :: statid
  character(len=16)  :: rsnumber
  character(len=16)  :: rssoftware

  call codes_open_file(ifile, '../../data/bufr/PraticaTemp.bufr', 'r')

  ! the first BUFR message is loaded from file
  ! ibufr is the BUFR id to be used in subsequent calls
  call codes_bufr_new_from_file(ifile, ibufr, iret)

  ! loop through all messages in the file
  do while (iret /= CODES_END_OF_FILE .AND. status_time == CODES_SUCCESS)

    ! Can check the template used
    ! call codes_get(ibufr,'unexpandedDescriptors',descriptors)
    ! write(0,*) 'Template: ', descriptors
    ! IF (descriptors(1) /= 309056.0) GOTO 999 ! only list descent profiles

    ! we need to instruct ecCodes to expand all the descriptors
    ! i.e. unpack the data values
    call codes_set(ibufr, "unpack", 1); 
    ! In our BUFR message verticalSoundingSignificance is always followed by
    !      geopotential, airTemperature, dewpointTemperature,
    !      windDirection, windSpeed and pressure.

    count = count + 1
    llskip = .False.

    ! Metadata:
    call codes_get(ibufr, 'shipOrMobileLandStationIdentifier', statid, status_id)
    IF (status_id /= CODES_SUCCESS) statid = "UNKNOWN"
    call codes_is_missing(ibufr, 'shipOrMobileLandStationIdentifier', statid_missing)
    IF (statid_missing == 1) statid = "MISSING"
    call codes_get(ibufr, 'blockNumber', blockNumber)
    call codes_get(ibufr, 'stationNumber', stationNumber)
    IF (blockNumber <= 99.0 .AND. stationNumber <= 1000) write (statid, '(I2.2,I3.3,3X)') blockNumber, stationNumber

    call codes_get(ibufr, 'year', year)
    call codes_get(ibufr, 'month', month)
    call codes_get(ibufr, 'day', day)
    call codes_get(ibufr, 'hour', hour)
    call codes_get(ibufr, 'minute', minute)
    call codes_get(ibufr, 'second', second, status_time)
    IF (status_time /= CODES_SUCCESS) second = 0.0
    call codes_get(ibufr, 'latitude', lat)
    call codes_get(ibufr, 'longitude', lon)
    call codes_get(ibufr, 'heightOfStationGroundAboveMeanSeaLevel', htg, status_ht)
    IF (status_ht /= CODES_SUCCESS) htg = -999.0
    call codes_get(ibufr, 'heightOfBarometerAboveMeanSeaLevel', htp, status_ht)
    IF (status_ht /= CODES_SUCCESS) htp = -999.0
    call codes_get(ibufr, 'radiosondeType', sondeType)
    call codes_get(ibufr, 'heightOfStation', htec, status_ht) ! Height from WMO list (BUFR)
    IF (status_ht == CODES_SUCCESS .AND. htg == -999.0) htg = htec
    ymd = INT(year)*10000 + INT(month)*100 + INT(day)
    hms = INT(hour)*10000 + INT(minute)*100 + INT(second)
    call codes_get(ibufr, 'radiosondeSerialNumber', rsnumber, status_rsno)
    call codes_get(ibufr, 'softwareVersionNumber', rssoftware, status_rssoft)
    call codes_get(ibufr, 'weightOfBalloon', balloonwt, status_balloonwt)
    IF (status_balloonwt /= CODES_SUCCESS) balloonwt = 0.0

    ! Ascent (skip incomplete reports for now)
    call codes_get(ibufr, 'timePeriod', timeVal, status_time)
    IF (status_time /= CODES_SUCCESS) THEN
      write (*, '(A,I7,A,A8,I9,I7.6,F9.3,F10.3,2F7.1,I4)') 'Ob: ', count, &
        ' ', statid, ymd, hms, lat(1), lon(1), htg, htp, INT(sondeType)
      write (*, '(A)') 'Missing times - skip'
      llskip = .True.
    END IF
    call codes_get(ibufr, 'pressure', presVal, status_p)
    IF (status_p /= CODES_SUCCESS) THEN
      write (*, '(A,I7,A,A8,I9,I7.6,F9.3,F10.3,2F7.1,I4)') 'Ob: ', count, &
        ' ', statid, ymd, hms, lat(1), lon(1), htg, htp, INT(sondeType)
      write (*, '(A)') 'Missing pressures - skip'
      llskip = .True.
    END IF
    call codes_get(ibufr, 'nonCoordinateGeopotentialHeight', zVal, status_ht)
    IF (status_ht /= CODES_SUCCESS) THEN
      write (*, '(A,I7,A,A8,I9,I7.6,F9.3,F10.3,2F7.1,I4)') 'Ob: ', count, &
        ' ', statid, ymd, hms, lat(1), lon(1), htg, htp, INT(sondeType)
      write (*, '(A)') 'Missing heights - skip'
      llskip = .True.
    END IF
    ! IF (blockNumber /= 17 .OR. stationNumber /= 196) llskip=.True.  ! FIX
    ! IF (blockNumber /= 17.0) llskip=.True.  ! FIX

    IF (.NOT. llskip) THEN
      call codes_get(ibufr, 'latitudeDisplacement', dlatVal)
      call codes_get(ibufr, 'longitudeDisplacement', dlonVal)
      call codes_get(ibufr, 'extendedVerticalSoundingSignificance', vssVal)
      call codes_get(ibufr, 'airTemperature', tVal)
      call codes_get(ibufr, 'dewpointTemperature', tdVal)
      call codes_get(ibufr, 'windDirection', wdirVal)
      call codes_get(ibufr, 'windSpeed', wspVal)

      ! ---- Array sizes (pressure size can be larger - wind shear levels)
      sizews = size(wspVal)

      ! ---- Print the values --------------------------------
      write (*, '(A,I7,A,A8,I9,I7.6,F9.3,F10.3,2F7.1,I4,I5)') 'Ob: ', count, &
        ' ', statid, ymd, hms, lat(1), lon(1), htg, htp, INT(sondeType), sizews
      IF (status_rsno == CODES_SUCCESS) write (*, '(A,A,A,F7.3)') &
        'RS number/software/balloonwt: ', rsnumber, rssoftware, balloonwt
      IF (status_ht == CODES_SUCCESS .AND. SIZE(lat) > 1) write (*, '(A,F9.3,F10.3,F7.1)') &
        'WMO list lat, lon, ht: ', lat(2), lon(2), htec
      write (*, '(A)') 'level  dtime   dlat   dlon pressure geopotH airTemp  dewPtT windDir  windSp  signif'
      do i = 1, sizews
        iflag = vssVal(i)
        IF (.NOT. llstdonly .OR. BTEST(iflag, 16)) &
          write (*, '(I5,F7.1,2F7.3,F9.1,F8.1,4F8.2,I8)') i, timeVal(i), &
          dlatVal(i), dlonVal(i), presVal(i), zVal(i), tVal(i), tdVal(i), &
          wdirVal(i), wspVal(i), INT(vssVal(i))
      end do

    END IF
    ! free allocated arrays
    IF (ALLOCATED(timeVal)) deallocate (timeVal)
    IF (ALLOCATED(dlatVal)) deallocate (dlatVal)
    IF (ALLOCATED(dlonVal)) deallocate (dlonVal)
    IF (ALLOCATED(vssVal)) deallocate (vssVal)
    IF (ALLOCATED(presVal)) deallocate (presVal)
    IF (ALLOCATED(zVal)) deallocate (zVal)
    IF (ALLOCATED(tVal)) deallocate (tVal)
    IF (ALLOCATED(tdVal)) deallocate (tdVal)
    IF (ALLOCATED(wdirVal)) deallocate (wdirVal)
    IF (ALLOCATED(wspVal)) deallocate (wspVal)
    IF (ALLOCATED(lat)) deallocate (lat)
    IF (ALLOCATED(lon)) deallocate (lon)

    ! 999 CONTINUE
    ! release the BUFR message
    call codes_release(ibufr)

    ! load the next BUFR message
    call codes_bufr_new_from_file(ifile, ibufr, iret)

  end do

  call codes_close_file(ifile)

end program bufr_read_tempf



Tabs Page
titleC


Code Block
languagecpp
titlebufr_read_tempf.c
/*
 * (C) Copyright 2005- ECMWF.
 *
 * This software is licensed under the terms of the Apache Licence Version 2.0
 * which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
 *
 * In applying this licence, ECMWF does not waive the privileges and immunities granted to it by
 * virtue of its status as an intergovernmental organisation nor does it submit to any jurisdiction.
 */

/*
 * C Implementation: bufr_read_tempf
 *
 * Description: read and print radiosonde data from TEMP BUFR messages.
 * If available this version also lists the position information from the WMO list
 * (now OSCAR/Surface) appended to the reports by ECMWF
 * 
 * Author: Bruce Ingleby
 */

/*
 * Please note that TEMP reports can be encoded in various ways in BUFR. Therefore the code
 * below might not work directly for other types of TEMP messages than the one used in the
 * example. It is advised to use bufr_dump to understand the structure of the messages.
 */

#include "eccodes.h"

/* Returns 1 if the bit at 'pos' in 'var' is set. The counting starts at 0 */
#define BTEST(var, pos) ((var) & (1 << (pos)))

/* Helper function to fill a double array with values of 'key'. Client must free memory */
static int get_double_array(codes_handle* h, const char* key, double** arr, size_t* size)
{
    CODES_CHECK(codes_get_size(h, key, size), 0);
    *arr = (double*)malloc(*size * sizeof(double));
    return codes_get_double_array(h, key, *arr, size);
}

/* Helper function to fill an integer (=long) array with values of 'key'. Client must free memory */
static int get_long_array(codes_handle* h, const char* key, long** arr, size_t* size)
{
    CODES_CHECK(codes_get_size(h, key, size), 0);
    *arr = (long*)malloc(*size * sizeof(long));
    return codes_get_long_array(h, key, *arr, size);
}

/* Reset dimension of input array to 'newsize' and fill with 'fillValue' */
static void realloc_and_fill(double** arr, size_t newsize, double fillValue)
{
    size_t i;
    free(*arr);
    *arr = (double*)malloc(newsize * sizeof(double));
    for(i=0; i<newsize; ++i) *arr[i] = fillValue;
}

int main(int argc, char* argv[])
{
    FILE* in = NULL;

    /* Message handle. Required in all the ecCodes calls acting on a message.*/
    codes_handle* h = NULL;

    double *lat = NULL, *lon = NULL, *presVal = NULL, *zVal = NULL;
    double *dlatVal = NULL, *dlonVal = NULL;
    double *tVal = NULL, *tdVal = NULL, *wspVal = NULL, *wdirVal = NULL;
    double htg, htp, htec = 0, balloonwt;
    int err         = 0;
    int status_rsno = 0, status_ht = 0, status_airt = 0, status_dewt = 0, status_p = 0;
    int count     = 0;
    int llskip    = 0;
    int llstdonly = 1; /* Set True to list standard levels only */
    size_t i = 0, size = 0, sizews = 0, sizelats = 0;
    long blockNumber, stationNumber;
    long year, month, day, hour, minute, second, ymd, hms, sondeType;
    long *timeVal = NULL, *vssVal = NULL;
    const char* infile = "../../data/bufr/PraticaTemp.bufr";
    char statid[128]   = {0,};
    char dropid[128]   = {0,};
    char rsnumber[16] = {0,};
    char rssoftware[16] = {0,};

    in = fopen(infile, "rb");
    if (!in) {
        fprintf(stderr, "Error: unable to open file %s\n", infile);
        return 1;
    }

    /* Loop over the messages in the BUFR file */
    while ((h = codes_handle_new_from_file(NULL, in, PRODUCT_BUFR, &err)) != NULL || err != CODES_SUCCESS) {
        count++;
        if (h == NULL) {
            fprintf(stderr, "Error: unable to create handle for message %d\n", count);
            continue;
        }

        /* We need to instruct ecCodes to expand the descriptors i.e., unpack the data values */
        CODES_CHECK(codes_set_long(h, "unpack", 1), 0);

        /* In our BUFR message verticalSoundingSignificance is always followed by
         *    geopotential, airTemperature, dewpointTemperature,
         *    windDirection, windSpeed and pressure */
        llskip = 0;

        size = 1024;
        err  = codes_get_string(h, "aircraftRegistrationNumberOrOtherIdentification", dropid, &size);
        if (err) strcpy(dropid, "UNKNOWN");
        size = 1024;
        err  = codes_get_string(h, "shipOrMobileLandStationIdentifier", statid, &size);
        if (err) strcpy(statid, dropid);

        CODES_CHECK(codes_get_long(h, "blockNumber", &blockNumber), 0);
        CODES_CHECK(codes_get_long(h, "stationNumber", &stationNumber), 0);
        if (blockNumber < 99 && stationNumber < 1000)
            sprintf(statid, "%ld%ld", blockNumber, stationNumber);
        CODES_CHECK(codes_get_long(h, "year", &year), 0);
        CODES_CHECK(codes_get_long(h, "month", &month), 0);
        CODES_CHECK(codes_get_long(h, "day", &day), 0);
        CODES_CHECK(codes_get_long(h, "hour", &hour), 0);
        CODES_CHECK(codes_get_long(h, "minute", &minute), 0);
        err = codes_get_long(h, "second", &second);
        if (err) second = 0;

        err = get_double_array(h, "latitude", &lat, &sizelats);

        err = get_double_array(h, "longitude", &lon, &size);

        status_ht = codes_get_double(h, "heightOfStationGroundAboveMeanSeaLevel", &htg);
        if (status_ht) htg = -999.0;
        status_ht = codes_get_double(h, "heightOfBarometerAboveMeanSeaLevel", &htp);
        if (status_ht) htp = -999.0;

        CODES_CHECK(codes_get_long(h, "radiosondeType", &sondeType), 0);
        err = codes_get_double(h, "heightOfStation", &htec); /* Height from WMO list (BUFR) */
        if (!err && htg == -999.0) htg = htec;

        ymd = year * 10000 + month * 100 + day;
        hms = hour * 10000 + minute * 100 + second;

        size        = 16;
        status_rsno = codes_get_string(h, "radiosondeSerialNumber", rsnumber, &size);
        size        = 16;
        codes_get_string(h, "softwareVersionNumber", rssoftware, &size);

        err = codes_get_double(h, "weightOfBalloon", &balloonwt);
        if (err) balloonwt = 0;

        /* Ascent (skip reports without dtime array for now) */
        err     = get_long_array(h, "timePeriod", &timeVal, &size);
        if (err) {
            printf("Ob: %d %s %ld %ld %g %g %g %g %ld\n",
                   count, statid, ymd, hms, lat[0], lon[0], htg, htp, sondeType);
            printf("Missing times - skip\n");
            llskip = 1;
        }
        status_p = get_double_array(h, "pressure", &presVal, &size);
        status_ht = get_double_array(h, "nonCoordinateGeopotentialHeight", &zVal, &size);

        if (!llskip) {
            err = get_double_array(h, "latitudeDisplacement", &dlatVal, &size);
            err = get_double_array(h, "longitudeDisplacement", &dlonVal, &size);
            err = get_long_array(h, "extendedVerticalSoundingSignificance", &vssVal, &size);
            status_airt = get_double_array(h, "airTemperature", &tVal, &size);
            status_dewt = get_double_array(h, "dewpointTemperature", &tdVal, &size);
            err = get_double_array(h, "windDirection", &wdirVal, &size);
            err = get_double_array(h, "windSpeed", &wspVal, &sizews);

            if (status_p != CODES_SUCCESS) {
                realloc_and_fill(&presVal, sizews, -999999999.0);
            }
            if (status_ht != CODES_SUCCESS) {
                realloc_and_fill(&zVal, sizews, -999999999.0);
            }
            if (status_airt != CODES_SUCCESS) {
                realloc_and_fill(&tVal, sizews, -999999999.0);
            }
            if (status_dewt != CODES_SUCCESS) {
                realloc_and_fill(&tdVal, sizews, -999999999.0);
            }
            /* Print the values */
            printf("Ob: %7d %s %ld %ld %7.3f %7.3f %7.1f %7.1f %4ld %5lu\n",
                   count, statid, ymd, hms, lat[0], lon[0], htg, htp, sondeType, sizews);
            if (status_rsno == CODES_SUCCESS) {
                printf("RS number/software/balloonwt: %s %s %7.3f\n", rsnumber, rssoftware, balloonwt);
            }
            if (status_ht == CODES_SUCCESS && sizelats > 1) {
                printf("WMO list lat, lon, ht: %s %g %g %g\n", statid, lat[1], lon[1], htec);
            }
            printf("level  dtime    dlat   dlon   pressure  geopotH  airTemp   dewPtT  windDir   windSp   signif\n");
            for (i = 0; i < sizews; ++i) {
                long iflag = vssVal[i];
                if (!llstdonly || BTEST(iflag, 16)) {
                    printf("%5lu %6ld %7.3f %7.3f %9.1f %8.1f %8.2f %8.2f %8.2f %8.2f %8ld\n",
                           i + 1, timeVal[i],
                           dlatVal[i], dlonVal[i],
                           presVal[i], zVal[i], tVal[i], tdVal[i],
                           wdirVal[i], wspVal[i], vssVal[i]);
                }
            }
        }

        /* Release memory */
        free(lat);
        free(lon);
        free(timeVal);
        free(dlatVal);
        free(dlonVal);
        free(presVal);
        free(zVal);
        free(tVal);
        free(tdVal);
        free(wdirVal);
        free(wspVal);
        free(vssVal);
        codes_handle_delete(h);
    }

    fclose(in);
    printf("Finishing normally. Number of BUFR records read: %d\n", count);
    return 0;
}




Tabs Page
titlePython


Code Block
languagepy
titlebufr_read_tempf.py
#
# (C) Copyright 2005- ECMWF.
#
# This software is licensed under the terms of the Apache Licence Version 2.0
# which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
#
# In applying this licence, ECMWF does not waive the privileges and immunities
# granted to it by virtue of its status as an intergovernmental organisation
# nor does it submit to any jurisdiction.
#
# Python implementation: bufr_read_tempf
#
#
# Description: read and print radiosonde data from TEMP BUFR messages.
# If available this version also lists the position information from the WMO list
# (now OSCAR/Surface) appended to the reports by ECMWF
#
# Author: Bruce Ingleby
#
# Please note that TEMP reports can be encoded in various ways in BUFR.
# Therefore the code below might not work directly for other types of TEMP
# messages than the one used in the example. It is advised to use bufr_dump to
# understand the structure of the messages.
#
from __future__ import print_functionsys
import traceback
import sys
import numpy as np

from eccodes import *

INPUT = "../../data/bufr/PraticaTemp.bufr"
VERBOSE = 1  # verbose error reporting


def example():
    # open BUFR file
    f = open(INPUT, "rb")
    llstdonly = 1
    cnt = 0
    # loop over the messages in the file
    while 1:
        # get handle for message
        bufr = codes_bufr_new_from_file(f)
        if bufr is None:
            break
        cnt += 1

        # desc = codes_get_array(bufr, 'unexpandedDescriptors')
        # if all(desc != 309056):    # descent reports
        #     codes_release(bufr)
        #     continue   # Skip other templates

        # we need to instruct ecCodes to expand all the descriptors
        # i.e. unpack the data section
        codes_set(bufr, "unpack", 1)
        # get header information from the message
        try:
            sid = codes_get(bufr, "aircraftRegistrationNumberOrOtherIdentification")
        except Exception:
            sid = "UNKNOWN"

        statid = "00000   "
        try:
            block = codes_get(bufr, "blockNumber")
            stnum = codes_get(bufr, "stationNumber")
            if (block > 0) and (block < 100):  # or block != CODES_MISSING_LONG
                statid = str.format("%.2i%.3i   " % (block, stnum))
        except Exception:
            statid = "00000   "
        if statid == "00000   ":
            statid = sid[0:8]

        # subtype = codes_get(bufr,'rdbSubtype')
        sondetype = codes_get(bufr, "radiosondeType")
        slat = codes_get_array(bufr, "latitude")
        slon = codes_get_array(bufr, "longitude")
        try:
            htg = codes_get(bufr, "heightOfStationGroundAboveMeanSeaLevel")
        except Exception:
            htg = -999.0
        try:
            htp = codes_get(bufr, "heightOfBarometerAboveMeanSeaLevel")
        except Exception:
            htp = -999.0
        year = codes_get(bufr, "year")
        month = codes_get(bufr, "month")
        day = codes_get(bufr, "day")
        hour = codes_get(bufr, "hour")
        minute = codes_get(bufr, "minute")
        try:
            second = codes_get(bufr, "second")
        except Exception:
            second = 0.0
        date = str.format("%i%.2i%.2i" % (year, month, day))
        time = str.format("%.2i%.2i%.2i" % (hour, minute, second))
        try:
            windsp = codes_get_array(bufr, "windSpeed")
        except Exception:
            codes_release(bufr)
            continue
        print(
            "Ob: %7i %s %s %s %7.3f %8.3f %6.1f %6.1f %3i %4i"
            % (
                cnt,
                statid,
                date,
                time,
                slat[0],
                slon[0],
                htg,
                htp,
                sondetype,
                len(windsp),
            )
        )

        try:
            rsnumber = codes_get(bufr, "radiosondeSerialNumber")
            rssoftware = codes_get(bufr, "softwareVersionNumber")
            balloonwt = codes_get(bufr, "weightOfBalloon")
            print("RS number/software/balloonwt ", rsnumber, rssoftware, balloonwt)
        except Exception:
            rsnumber = 0
        try:
            htec = codes_get(
                bufr, "heightOfStation"
            )  # Height from WMO list (appended by ECMWF)
            print("WMO list lat, lon, ht: %7.3f %8.3f %6.1f" % (slat[1], slon[1], htec))
        except Exception:
            htec = 0

        # get all the timePeriods
        dtime = codes_get_array(bufr, "timePeriod")
        try:
            pressure = codes_get_array(bufr, "pressure")
        except Exception:
            codes_release(bufr)
            continue
        vsSignif = codes_get_array(bufr, "extendedVerticalSoundingSignificance")
        try:
            geopoth = codes_get_array(bufr, "nonCoordinateGeopotentialHeight")
        except Exception:
            codes_release(bufr)
            continue
        dlat = codes_get_array(bufr, "latitudeDisplacement")
        dlon = codes_get_array(bufr, "longitudeDisplacement")
        airt = codes_get_array(bufr, "airTemperature")
        dewt = codes_get_array(bufr, "dewpointTemperature")
        windd = codes_get_array(bufr, "windDirection")
        dtime = np.where(dtime != CODES_MISSING_LONG, dtime, np.nan)
        dlat = np.where(dlat != CODES_MISSING_DOUBLE, dlat, np.nan)
        dlon = np.where(dlon != CODES_MISSING_DOUBLE, dlon, np.nan)
        airt = np.where(airt != CODES_MISSING_DOUBLE, airt, np.nan)
        dewt = np.where(dewt != CODES_MISSING_DOUBLE, dewt, np.nan)
        windd = np.where(windd != CODES_MISSING_LONG, windd, np.nan)
        windsp = np.where(windsp != CODES_MISSING_DOUBLE, windsp, np.nan)
        geopoth = np.where(geopoth != CODES_MISSING_DOUBLE, geopoth, np.nan)
        pressure = np.where(pressure != CODES_MISSING_DOUBLE, pressure, np.nan)
        # pressure = np.where(pressure > -1e10, pressure, np.nan)
        print(
            "level  dtime   dlat   dlon pressure geopotH airTemp  dewPtT windDir  windSp  signif"
        )
        for i in range(0, len(windsp)):
            if (not llstdonly) or vsSignif[i] != 65536:
                continue
            print(
                "%5i %6.1f %6.3f %6.3f %8.1f %7.1f %7.2f %7.2f %7.2f %7.2f %7i"
                % (
                    i + 1,
                    dtime[i],
                    dlat[i],
                    dlon[i],
                    pressure[i],
                    geopoth[i],
                    airt[i],
                    dewt[i],
                    windd[i],
                    windsp[i],
                    vsSignif[i],
                )
            )
        # delete handle
        codes_release(bufr)
    # close the file
    f.close()


def main():
    try:
        example()
    except CodesInternalError as err:
        if VERBOSE:
            traceback.print_exc(file=sys.stderr)
        else:
            sys.stderr.write(err.msg + "\n")
        return 1


if __name__ == "__main__":
    sys.exit(main())




...