NCEPLIBS-g2  3.4.7
gribend.F90
Go to the documentation of this file.
1 
5 
26 subroutine gribend(cgrib, lcgrib, lengrib, ierr)
27  implicit none
28 
29  character(len = 1), intent(inout) :: cgrib(lcgrib)
30  integer, intent(in) :: lcgrib
31  integer, intent(out) :: lengrib, ierr
32 
33  integer ilen, isecnum
34  character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
35  character(len = 4):: ctemp
36  integer iofst, lencurr, len
37 
38  ierr = 0
39 
40  ! Check to see if beginning of GRIB message exists.
41  ctemp = cgrib(1) // cgrib(2) // cgrib(3) // cgrib(4)
42  if (ctemp .ne. grib) then
43  print *, 'gribend: GRIB not found in given message.'
44  ierr = 1
45  return
46  endif
47 
48  ! Get current length of GRIB message.
49  call g2_gbytec(cgrib, lencurr, 96, 32)
50 
51  ! Loop through all current sections of the GRIB message to
52  ! find the last section number.
53  len = 16 ! Length of Section 0
54  do
55  ! Get number and length of next section.
56  iofst = len * 8
57  call g2_gbytec(cgrib, ilen, iofst, 32)
58  iofst = iofst + 32
59  call g2_gbytec(cgrib, isecnum, iofst, 8)
60  len = len + ilen
61 
62  ! Exit loop if last section reached.
63  if (len .eq. lencurr) exit
64 
65  ! If byte count for each section doesn't match current total
66  ! length, then there is a problem.
67  if (len .gt. lencurr) then
68  print *, 'gribend: Section byte counts don''t add ' &
69  ,'to total.'
70  print *, 'gribend: Sum of section byte counts = ', len
71  print *, 'gribend: Total byte count in Section 0 = ', &
72  lencurr
73  ierr = 3
74  return
75  endif
76  enddo
77 
78  ! Can only add End Section (Section 8) after Section 7.
79  if (isecnum .ne. 7) then
80  print *, 'gribend: Section 8 can only be added after Section 7.'
81  print *, 'gribend: Section ', isecnum, &
82  ' was the last found in',' given GRIB message.'
83  ierr = 4
84  return
85  endif
86 
87  ! Add Section 8 - End Section.
88  cgrib(lencurr + 1:lencurr + 4) = c7777
89 
90  ! Update current byte total of message in Section 0.
91  lengrib = lencurr + 4
92  call g2_sbytec(cgrib, lengrib, 96, 32)
93 end subroutine gribend
subroutine g2_gbytec(in, iout, iskip, nbits)
Extract arbitrary size values from a packed bit string, right justifying each value in the unpacked a...
Definition: g2_gbytesc.F90:20
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 gribend(cgrib, lcgrib, lengrib, ierr)
Finalize a GRIB2 message after all grids and fields have been added.
Definition: gribend.F90:27