NCEPLIBS-g2  3.4.5
addlocal.f
Go to the documentation of this file.
1 
6 
30 
31  subroutine addlocal(cgrib,lcgrib,csec2,lcsec2,ierr)
32 
33  character(len=1),intent(inout) :: cgrib(lcgrib)
34  character(len=1),intent(in) :: csec2(lcsec2)
35  integer,intent(in) :: lcgrib,lcsec2
36  integer,intent(out) :: ierr
37 
38  character(len=4),parameter :: grib='GRIB',c7777='7777'
39  character(len=4):: ctemp
40  integer,parameter :: two=2
41  integer lensec2,iofst,ibeg,lencurr,len
42 
43  ierr=0
44 
45 ! Check to see if beginning of GRIB message exists
46 
47  ctemp=cgrib(1)//cgrib(2)//cgrib(3)//cgrib(4)
48  if ( ctemp.ne.grib ) then
49  print *,'addlocal: GRIB not found in given message.'
50  print *,'addlocal: Call to routine gribcreate required',
51  & ' to initialize GRIB messge.'
52  ierr=1
53  return
54  endif
55 
56 ! Get current length of GRIB message
57 
58  call g2_gbytec(cgrib,lencurr,96,32)
59 
60 ! Check to see if GRIB message is already complete
61 
62  ctemp=cgrib(lencurr-3)//cgrib(lencurr-2)//cgrib(lencurr-1)
63  & //cgrib(lencurr)
64  if ( ctemp.eq.c7777 ) then
65  print *,'addlocal: GRIB message already complete. Cannot',
66  & ' add new section.'
67  ierr=2
68  return
69  endif
70 
71 ! Loop through all current sections of the GRIB message to
72 ! find the last section number.
73 
74  len=16 ! length of Section 0
75  do
76 ! Get section number and length of next section
77  iofst=len*8
78  call g2_gbytec(cgrib,ilen,iofst,32)
79  iofst=iofst+32
80  call g2_gbytec(cgrib,isecnum,iofst,8)
81  len=len+ilen
82 ! Exit loop if last section reached
83  if ( len.eq.lencurr ) exit
84 ! If byte count for each section doesn't match current
85 ! total length, then there is a problem.
86  if ( len.gt.lencurr ) then
87  print *,'addlocal: Section byte counts don''t add to total.'
88  print *,'addlocal: Sum of section byte counts = ',len
89  print *,'addlocal: Total byte count in Section 0 = ',lencurr
90  ierr=3
91  return
92  endif
93  enddo
94 
95 ! Section 2 can only be added after sections 1 and 7.
96 
97  if ( (isecnum.ne.1) .and. (isecnum.ne.7) ) then
98  print *,'addlocal: Section 2 can only be added after Section',
99  & ' 1 or Section 7.'
100  print *,'addlocal: Section ',isecnum,' was the last found in',
101  & ' given GRIB message.'
102  ierr=4
103  return
104  endif
105 
106 ! Add Section 2 - Local Use Section
107 
108  ibeg=lencurr*8 ! Calculate offset for beginning of section 2
109  iofst=ibeg+32 ! leave space for length of section
110  call g2_sbytec(cgrib,two,iofst,8) ! Store section number ( 2 )
111  istart=lencurr+5
112  cgrib(istart+1:istart+lcsec2)=csec2(1:lcsec2)
113 
114 ! Calculate length of section 2 and store it in octets
115 ! 1-4 of section 2.
116 
117  lensec2=lcsec2+5 ! bytes
118  call g2_sbytec(cgrib,lensec2,ibeg,32)
119 
120 
121 ! Update current byte total of message in Section 0
122 
123  call g2_sbytec(cgrib,lencurr+lensec2,96,32)
124 
125  return
126  end
127 
128 
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
addlocal
subroutine addlocal(cgrib, lcgrib, csec2, lcsec2, ierr)
This subroutine adds a Local Use Section (Section 2) to a GRIB2 message.
Definition: addlocal.f:32
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