Description

This example shows: how to set a bitmap in a GRIB message to encode missing values in the data

Source code

grib_set_bitmap.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: how to set a bitmap in a GRIB message
!               to encode missing values in the data
!
!
program set_bitmap
  use eccodes
  implicit none
  integer                         :: infile,outfile
  integer                         :: igrib, iret
  integer                         :: numberOfValues
  real, dimension(:), allocatable :: values
  real                            :: missingValue
  logical                         :: grib1Example

  grib1Example=.true.

  if (grib1Example) then
    ! GRIB 1 example
    call codes_open_file(infile,'../../data/regular_latlon_surface.grib1','r')
  else
    ! GRIB 2 example
    call codes_open_file(infile,'../../data/regular_latlon_surface.grib2','r')
  end if

  call codes_open_file(outfile,'out.set_bitmap_f.grib','w')

  ! A new grib message is loaded from file
  ! igrib is the grib id to be used in subsequent calls
  call codes_grib_new_from_file(infile,igrib)

  ! The missingValue is not coded in the message.
  ! It is a value we define as a placeholder for a missing value
  ! at a point in the grid.
  ! It should be chosen so that it cannot be confused
  ! with a valid field value
  missingValue=1.0e36
  call codes_set(igrib, 'missingValue',missingValue)
  write(*,*) 'missingValue=',missingValue

  ! get the size of the values array
  call codes_get_size(igrib,'values',numberOfValues)
  write(*,*) 'numberOfValues=',numberOfValues

  allocate(values(numberOfValues), stat=iret)

  ! get data values
  call codes_get(igrib,'values',values)

  ! enable bitmap
  call codes_set(igrib, 'bitmapPresent', 1)

  ! set some values to be missing
  values(1:10) = missingValue

  ! set the values (the bitmap will be automatically built)
  call codes_set(igrib,'values', values)

  !  write modified message to a file
  call codes_write(igrib,outfile)

  ! free memory
  call codes_release(igrib)

  call codes_close_file(infile)
  call codes_close_file(outfile)

  deallocate(values)

end program set_bitmap
grib_set_bitmap.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.
#
# Description: how to set a bitmap in a GRIB message
#              to encode missing values in the data
#
import traceback
import sys
from eccodes import *

INPUT = '../../data/regular_latlon_surface.grib1'
OUTPUT = 'out.set_bitmap_p.grib'
MISSING = 1.0e36
VERBOSE = 1  # verbose error reporting


def example():
    fin = open(INPUT, 'rb')
    fout = open(OUTPUT, 'wb')
    gid = codes_grib_new_from_file(fin)

    # The missingValue is not coded in the message.
    # It is a value we define as a placeholder for a missing value
    # at a point in the grid.
    # It should be chosen so that it cannot be confused
    # with a valid field value
    codes_set(gid, 'missingValue', MISSING)

    values = codes_get_values(gid)

    # Enable bitmap
    codes_set(gid, 'bitmapPresent', 1)

    # Change some data values to be missing
    num_missing = 0
    for i in range(100):
        # Set every other value to a missing value
        if i % 2 == 0:
            values[i] = MISSING
            num_missing += 1
    codes_set_values(gid, values)

    # Check counts of missing and non-missing values
    num_data = codes_get(gid, 'numberOfDataPoints', int)
    assert num_data == len(values)
    assert codes_get(gid, 'numberOfCodedValues', int) == num_data - num_missing
    assert codes_get(gid, 'numberOfMissing', int) == num_missing

    codes_write(gid, fout)
    codes_release(gid)
    fin.close()
    fout.close()


def main():
    try:
        example()
    except CodesInternalError as err:
        if VERBOSE:
            traceback.print_exc(file=sys.stderr)
        else:
            print(err.msg, file=sys.stderr)

        return 1


if __name__ == "__main__":
    sys.exit(main())
grib_set_bitmap.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: grib_set_bitmap
 *
 * Description: how to set a bitmap in a GRIB message
 *              to encode missing values in the data
 *
 */

#include <stdio.h>
#include <stdlib.h>

#include "eccodes.h"

int main(int argc, char** argv)
{
    int err = 0;
    size_t size=0;

    FILE* in = NULL;
    const char* infile = "../../data/regular_latlon_surface.grib1";
    FILE* out = NULL;
    const char* outfile = "out.set_bitmap_c.grib";
    codes_handle *h = NULL;
    const void* buffer = NULL;
    size_t values_len;
    double* values;
    double missing=1.0e36;
    int i=0;

    in = fopen(infile, "r");
    if(!in) {
        printf("ERROR: unable to open input file %s\n",infile);
        return 1;
    }

    out = fopen(outfile, "w");
    if(!out) {
        printf("ERROR: unable to open output file %s\n",outfile);
        fclose(in);
        return 1;
    }

    h = codes_handle_new_from_file(0, in, PRODUCT_GRIB, &err);
    if (h == NULL) {
        printf("Error: unable to create handle from file %s\n",infile);
    }

    CODES_CHECK(codes_set_double(h, "missingValue", missing),0);

    /* get the size of the values array*/
    CODES_CHECK(codes_get_size(h, "values", &values_len),0);

    values = (double*)malloc(values_len*sizeof(double));

    /* get data values*/
    CODES_CHECK(codes_get_double_array(h, "values", values, &values_len),0);

    CODES_CHECK(codes_set_long(h, "bitmapPresent", 1),0);

    for(i = 0; i < 10; i++) {
        values[i]=missing;
    }

    CODES_CHECK(codes_set_double_array(h, "values", values, values_len),0);

    /* get the coded message in a buffer */
    CODES_CHECK(codes_get_message(h, &buffer, &size),0);

    /* write the buffer in a file*/
    if(fwrite(buffer, 1, size, out) != size) {
        perror(outfile);
        exit(1);
    }

    codes_handle_delete(h);

    fclose(in);
    fclose(out);

    return 0;
}