NCEPLIBS-g2  3.4.5
gribend.f
Go to the documentation of this file.
1 
6 
30 
31  subroutine gribend(cgrib,lcgrib,lengrib,ierr)
32 
33  character(len=1),intent(inout) :: cgrib(lcgrib)
34  integer,intent(in) :: lcgrib
35  integer,intent(out) :: lengrib,ierr
36 
37  character(len=4),parameter :: grib='GRIB',c7777='7777'
38  character(len=4):: ctemp
39  integer iofst,ibeg,lencurr,len
40 
41  ierr=0
42 !
43 ! Check to see if beginning of GRIB message exists
44 !
45  ctemp=cgrib(1)//cgrib(2)//cgrib(3)//cgrib(4)
46  if ( ctemp.ne.grib ) then
47  print *,'gribend: GRIB not found in given message.'
48  ierr=1
49  return
50  endif
51 !
52 ! Get current length of GRIB message
53 !
54  call g2_gbytec(cgrib,lencurr,96,32)
55 !
56 ! Check to see if GRIB message is already complete
57 !
58 ! ctemp=cgrib(lencurr-3)//cgrib(lencurr-2)//cgrib(lencurr-1)
59 ! & //cgrib(lencurr)
60 ! if ( ctemp.eq.c7777 ) then
61 ! print *,'gribend: GRIB message already complete.'
62 ! ierr=2
63 ! return
64 ! endif
65 !
66 ! Loop through all current sections of the GRIB message to
67 ! find the last section number.
68 !
69  len=16 ! Length of Section 0
70  do
71  ! Get number and length of next section
72  iofst=len*8
73  call g2_gbytec(cgrib,ilen,iofst,32)
74  iofst=iofst+32
75  call g2_gbytec(cgrib,isecnum,iofst,8)
76  len=len+ilen
77  ! Exit loop if last section reached
78  if ( len.eq.lencurr ) exit
79  ! If byte count for each section doesn't match current
80  ! total length, then there is a problem.
81  if ( len.gt.lencurr ) then
82  print *,'gribend: Section byte counts don''t add to total.'
83  print *,'gribend: Sum of section byte counts = ',len
84  print *,'gribend: Total byte count in Section 0 = ',lencurr
85  ierr=3
86  return
87  endif
88  enddo
89 !
90 ! Can only add End Section (Section 8) after Section 7.
91 !
92  if ( isecnum.ne.7 ) then
93  print *,'gribend: Section 8 can only be added after Section 7.'
94  print *,'gribend: Section ',isecnum,' was the last found in',
95  & ' given GRIB message.'
96  ierr=4
97  return
98  endif
99 !
100 ! Add Section 8 - End Section
101 !
102  cgrib(lencurr+1:lencurr+4)=c7777
103 
104 !
105 ! Update current byte total of message in Section 0
106 !
107  lengrib=lencurr+4
108  call g2_sbytec(cgrib,lengrib,96,32)
109 
110  return
111  end
112 
113 
114 
115 
g2_sbytec
subroutine g2_sbytec(OUT, IN, ISKIP, NBYTE)
This subrountine is to put arbitrary size values into a packed bit string, taking the low order bits ...
Definition: g2_gbytesc.f:39
g2_gbytec
subroutine g2_gbytec(IN, IOUT, ISKIP, NBYTE)
This subrountine is to extract arbitrary size values from a packed bit string, right justifying each ...
Definition: g2_gbytesc.f:20
gribend
subroutine gribend(cgrib, lcgrib, lengrib, ierr)
This subroutine finalizes a GRIB message after all grids and fields have been added.
Definition: gribend.f:32