NCEPLIBS-g2  3.4.7
gribcreate.F90
Go to the documentation of this file.
1 
5 
54 subroutine gribcreate(cgrib, lcgrib, listsec0, listsec1, ierr)
55  implicit none
56 
57  character(len = 1), intent(inout) :: cgrib(lcgrib)
58  integer, intent(in) :: listsec0(*), listsec1(*)
59  integer, intent(in) :: lcgrib
60  integer, intent(out) :: ierr
61 
62  integer :: i, lensec1, nbits
63  character(len = 4), parameter :: grib = 'GRIB'
64  integer, parameter :: ZERO = 0, one = 1
65  integer, parameter :: mapsec1len = 13
66  integer, parameter :: mapsec1(mapsec1len) = (/ 2, 2, 1, 1, 1, 2, 1, 1, 1, 1, 1, 1, 1 /)
67  integer lensec0, iofst, ibeg
68 
69  ierr = 0
70 
71 #ifdef LOGGING
72  print *, 'gribcreate lcgrib ', lcgrib
73 #endif
74 
75  ! Currently handles only GRIB Edition 2.
76  if (listsec0(2) .ne. 2) then
77  print *, 'gribcreate: can only code GRIB edition 2.'
78  ierr = 1
79  return
80  endif
81 
82  ! Pack Section 0 - Indicator Section (except for total length of
83  ! GRIB message).
84  cgrib(1) = grib(1:1) ! Beginning of GRIB message
85  cgrib(2) = grib(2:2)
86  cgrib(3) = grib(3:3)
87  cgrib(4) = grib(4:4)
88  call g2_sbytec(cgrib, zero, 32, 16) ! reserved for future use
89  call g2_sbytec(cgrib, listsec0(1), 48, 8) ! Discipline
90  call g2_sbytec(cgrib, listsec0(2), 56, 8) ! GRIB edition number
91  lensec0 = 16 ! bytes (octets)
92 
93  ! Pack Section 1 - Identification Section.
94  ibeg = lensec0 * 8 ! Calculate offset for beginning of section 1
95  iofst = ibeg + 32 ! leave space for length of section
96  call g2_sbytec(cgrib, one, iofst, 8) ! Store section number ( 1 )
97  iofst = iofst + 8
98 
99  ! Pack up each input value in array listsec1 into the the
100  ! appropriate number of octets, which are specified in corresponding
101  ! entries in array mapsec1.
102  do i = 1, mapsec1len
103  nbits = mapsec1(i) * 8
104  call g2_sbytec(cgrib, listsec1(i), iofst, nbits)
105  iofst = iofst + nbits
106  enddo
107 
108  ! Calculate length of section 1 and store it in octets 1-4 of
109  ! section 1.
110  lensec1 = (iofst - ibeg) / 8
111  call g2_sbytec(cgrib, lensec1, ibeg, 32)
112 
113  ! Put current byte total of message into Section 0.
114  call g2_sbytec(cgrib, zero, 64, 32)
115  call g2_sbytec(cgrib, lensec0 + lensec1, 96, 32)
116 end subroutine gribcreate
subroutine g2_sbytec(out, in, iskip, nbits)
Put arbitrary size values into a packed bit string, taking the low order bits from each value in the ...
Definition: g2_gbytesc.F90:42
subroutine gribcreate(cgrib, lcgrib, listsec0, listsec1, ierr)
This subroutine initializes a new GRIB2 message and packs GRIB2 sections 0 (Indicator Section) and 1 ...
Definition: gribcreate.F90:55