This shell script uses a Fortran program calling ecCodes routines to decode GRIB data.

  • fields are retrieved from MARS in GRIB format
  • the GRIB decoding Fortran program is compiled and linked  against the ecCodes library using gfortran
  • the program is executed to decode and prints various sections of GRIB data using calls to the grib_api ecCodes functions


#!/bin/bash
#
# **************************** LICENSE START ***********************************
#
# Copyright 2021 ECMWF. This software is distributed under the terms
# of the Apache License version 2.0. In applying this license, ECMWF does not
# waive the privileges and immunities granted to it by virtue of its status as
# an Intergovernmental Organization or submit itself to any jurisdiction.
#
# ***************************** LICENSE END ************************************
#
# Retrieve_decode_grib_api      USER SERVICES  NOVEMBER 2021 - ECMWF
#
# 
#       This shell-script:
#
#          - retrieves fields in GRIB format from MARS
#          - compiles and links a Fortran90 decoding program using gfortran 
#          - decodes and prints various sections of GRIB data, using ecCodes
#
#       This shell script produces the standard output file 
#
#         Retrieve_decode_grib.<JOB-ID>.out
#
#       containing the log of job execution.    
#
#       For more information on ecCodes see 
#
#       https://confluence.ecmwf.int/display/ECC/ecCodes+Home
#
#
#-------------------------------
# setting options for SLURM
#-------------------------------
# Options that are specified within the script file should precede the
# first executable shell command in the file. 
# All options, which are case sensitive, are preceded by #SBATCH. 
# These lines become active only when the script is submitted
# using the "sbatch" command. 
# All job output is written to the workdir directory, by default. 


#SBATCH --qos=ef

        # Specifies that your job will run in the queue (Quality Of
        # Service) "ef".

#SBATCH --job-name=retrieve_decode_grib_api

        # Assigns the specified name to the request

#SBATCH --output=retrieve_decode_grib_api.%j.out

        # Specifies the name and location of STDOUT where %j is the job-id
        # The file will be # written in the workdir directory if it is a 
        # relative path. If not given, the default is slurm-%j.out in the
        # workdir.

#SBATCH --error=retrieve_decode_grib_api.%j.out

        # Specifies the name and location of STDERR where %j is the job-id
        # The file will be # written in the workdir directory if it is a 
        # relative path. If not given, the default is slurm-%j.outin the
        # workdir.

#SBATCH --chdir=/scratch/...

        # Sets the working directory of the batch script before it is
        # executed.

#SBATCH --mail-type=FAIL
        
        # Specifies that an email should be sent in case the job fails.
        # Other options include BEGIN, END, REQUEUE and ALL (any state
        # change).

#SBATCH --time=00:05:00

        # Specifies that your job my run up to HH:MM:SS of wall clock
        # time. The job will be killed if it exceeds this limit. If
        # time is not defined, the default limit for the queue (qos)
        # will be used.

#-------------------------------
# setting environment variables
#-------------------------------

export PATH=$PATH:.             # Allows you to run any of your programs or
                                # scripts held in the current directory (not 
                                # required if already done in your .user_profile 
                                # or .user_kshrc)
set -ev

#-------------------------------
# commands to be executed
#-------------------------------

cd $SCRATCHDIR		# All the files created in this directory will be
                        # deleted when the job terminates.


# cat is used to read the Fortran 90 program from the input stream 
# and to write it to the file 'grib_api_demo.f90'.
# cat reads line by line until it reaches a line which starts with EOF.
#

cat>grib_api_demo.f90<<EOF
PROGRAM grib_api_demo
  use grib_api
  implicit none

!                 

  INTEGER                          :: status, iret, rfile, gribid, i         
  CHARACTER(LEN=*), PARAMETER      :: input_file='grib_file.grb'
  CHARACTER(LEN=*), PARAMETER      :: open_mode='r'

  real, dimension(:), allocatable  ::  values
  integer                          ::  numberOfPoints


!                                                                       
!     Open data file for reading.                                        
!
  call grib_open_file (rfile, input_file, open_mode) 
!
!     Check return code.
!   
  call grib_new_from_file(rfile, gribid, status)

  LOOP: DO WHILE (status /= GRIB_END_OF_FILE)   

!**************************************************
!             Start of User data processing
!**************************************************

!
!   grib_dump should only be usd for diagnostic purposes.
!
!   To process the data and grib headers, you will have to use
!   grib_get and request the decoding of the keys you need,
!   e.g. the date, the time, the parameter, the "level".
!

    call grib_dump(gribid)
    call grib_get_size(gribid,'values',numberOfPoints)
    allocate(values(numberOfPoints), stat=iret)
    if (iret /= 0) then
       STOP 'grib_api_demo: allocate failed'
    end if
    call grib_get(gribid,'values',values)
    print *, "First 20 values:"
    write(*,'(f10.4)')(values(i),i=1,20)


!**************************************************
!             End of User data processing
!**************************************************
!
! To access the next record in the grib file
!
    call grib_release (gribid)

    deallocate(values)

    call grib_new_from_file(rfile, gribid, status)

!                 
!        Loop back for next input field.                                   
!            
        END DO LOOP
!                                                                          
!     End-of-file on input.             
!
! release and close grib file                     
  call grib_close_file (rfile)
!                                                             
!
END PROGRAM grib_api_demo    
EOF



#--------------------
# MARS  request
#--------------------

# Retrieve the data from MARS to the target data file "grib_file".

mars<<EOF
retrieve,
    class  = od,
    stream = oper,
    date   = -2,
    time   = 00,
    type   = fc,
    step   = 12,
    target = "grib_file.grb"
EOF


if [ $? != 0 ]                          # Check MARS exit code. 
then
      echo " The MARS request failed."
      exit 1
fi


#-------------------------------
# compile and link libraries
#-------------------------------

module load prgenv/gnu
module load ecmwf-toolbox

gfortran -march=native -O3 -o grib_api_demo grib_api_demo.f90 $ECCODES_INCLUDE $ECCODES_LIB       

                         # Use gfortran -o to create the grib_api_demo
                         # executable. The compiler options
                         # -arch=native -O3 provide
                         # some basic optimisation. The variables
                         # $ECCODES_INCLUDE and $ECCODES_LIB
                         # respectively contain the location of the
                         # ecCodes include files and libraries.

time ./grib_api_demo     # timing the run of grib_api_demo 


#-------------------------------
# tidy up by deleting unwanted files
#-------------------------------
# This is done automatically when using $SCRATCHDIR.

exit 0

# End of example job 'Retrieve_decode_grib_api'