Skip Navigation Links www.nws.noaa.gov 
NOAA logo - Click to go to the NOAA home page National Weather Service   NWS logo - Click to go to the NWS home page
Climate Prediction Center
 
 

 
About Us
   Our Mission
   Who We Are

Contact Us
   CPC Information
   CPC Web Team

 
HOME > Monitoring_and_Data > Oceanic and Atmospheric Data > Reanalysis: Atmospheric Data > wgrib2api: memory buffers/files
 

wgrib2api: memory buffers/files

Introduction

Wgrib2api can be used to read from and write to grib files which is great for non-HPC applications. For HPC, you want to parallelize the encoding and decoding of the grib messages. This can be done using memory buffers (memory files in wgrib2 terminology) instead of disk files. As a practical matter, wgrib2api does not support filelocks and only one instance of wgrib2api should be allowed to write to any disk file. Instead you have to use memory files where each instance of wgrib2api has its own independent copy of the memory files. For example, writing the inventory to memory file 0 can be done by,

   i = grb2_mk_inv('MYDATA.GRB', '@mem:0')         ! write index/inventory memory file 0
Suppose that you want to copy grib data to memory file 0, create an inventory using memory file 1 and then read the Z500 field from memory file 0. The code segment will go like,
    use wgrib2api
    real, allocatable :: grid(:,:)

!   character (len=1) buff(1000000)
!   buff contains a grib message

!   copy buff to memory buffer/file 0
    i = wgrib2_set_mem_buffer(buf, 1000000, 0)
    if (i.ne.0) stop 1

    i = grb2_mk_inv('@mem:0', '@mem:1')
    if (i.ne.0) stop 2

    i = grb2_read('@mem:0', '@mem:1', ':HGT:500 mb:', data2=grid)
    if (i.ne.1) stop 3

    write(*,*) ' read Z500 from buf(*) z(20,20)=', grid(20,20)
    stop
    end
You can encode grib to a memory file by,
    use wgrib2api
    real, allocatable :: grid(:,:)
    character, (len=1), allocatable :: buffer(:)

    nx=360
    ny=181
    allocate(grid(nx,ny))
    read(11) grid
    i = grb2_wrt('@mem:0','TEMPPLATE.grb',1,data2=grid,meta='D=20170102030000:HGT:500 mb:anl:')
    if (i .ne. 0) stop 1
    n = wgrib2_get_mem_buffer_size(0)
    allocate (buffer(n))
    i = wgrib2_get_mem_buffer(buffer, n, 0)

Usage

   SIZE = wgrib2_get_mem_buffer_size(FILE_NUM)
       returns the size of memory file '@mem:FILE_NUM'
       FILE_NUM:   integer, memory file number, 0..19
       SIZE:       integer, size of memory buffer

   i = wgrib2_get_mem_buffer(BUFFER, SIZE, FILE_NUM)
       copies memory file '@mem:FILE_NUM' to BUFFER(1:SIZE)
       BUFFER:     character (len=1) BUFFER(1:SIZE), read from memory file
       SIZE:       integer size to transfer, use get_mem_buffer_size(FILE_NUM)
       FILE_NUM:   integer, memory file number
       i:          integer, return code, 0 = worked

   i = wgrib2_set_mem_buffer(BUFFER, SIZE, FILE_NUM)
       copies BUFFER to memory file '@mem:FILE_NUM'
       BUFFER:     character (len=1) BUFFER(1:SIZE), write to memory file
       SIZE:       integer size to transfer
       FILE_NUM:   integer, memory file number
       i:          integer, return code, 0 = worked

Wgrib2api is made of a high level (grb2_*) functions and low_level functions. The low level routines can be found in wgrib2lowapi and are named wgrib2_*.


NOAA/ National Weather Service
National Centers for Environmental Prediction
Climate Prediction Center
5830 University Research Court
College Park, Maryland 20740
Climate Prediction Center Web Team
Page last modified: March 13, 2017
Disclaimer Privacy Policy