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