NCEPLIBS-g2  3.4.7
gettemplates.F90
Go to the documentation of this file.
1 
5 
66 subroutine gettemplates(cgrib, lcgrib, ifldnum, igds, igdstmpl, &
67  igdslen, ideflist, idefnum, ipdsnum, ipdstmpl, ipdslen, coordlist, numcoord, ierr)
68  implicit none
69 
70  character(len = 1), intent(in) :: cgrib(lcgrib)
71  integer, intent(in) :: lcgrib, ifldnum
72  integer, intent(out) :: igds(*), igdstmpl(*), ideflist(*)
73  integer, intent(out) :: ipdsnum, ipdstmpl(*)
74  integer, intent(out) :: idefnum, numcoord
75  integer, intent(out) :: ierr
76  real, intent(out) :: coordlist(*)
77 
78  character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
79  character(len = 4) :: ctemp
80  integer:: listsec0(2)
81  integer iofst, istart
82  logical have3, have4
83  integer :: igdslen, ipdslen, ipos, isecnum, j, jerr, lengrib, lensec, lensec0, numfld
84 
85  have3 = .false.
86  have4 = .false.
87  ierr = 0
88  numfld = 0
89 
90  ! Check for valid request number.
91  if (ifldnum .le. 0) then
92  print *, 'gettemplates: Request for field number must be positive.'
93  ierr = 3
94  return
95  endif
96 
97  ! Check for beginning of GRIB message in the first 100 bytes.
98  istart = 0
99  do j = 1, 100
100  ctemp = cgrib(j) // cgrib(j + 1) // cgrib(j + 2) // cgrib(j + 3)
101  if (ctemp .eq. grib ) then
102  istart = j
103  exit
104  endif
105  enddo
106  if (istart .eq. 0) then
107  print *, 'gettemplates: Beginning characters GRIB not found.'
108  ierr = 1
109  return
110  endif
111 
112  ! Unpack Section 0 - Indicator Section.
113  iofst = 8 * (istart + 5)
114  call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
115  iofst = iofst + 8
116  call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
117  iofst = iofst + 8
118  iofst = iofst + 32
119  call g2_gbytec(cgrib, lengrib, iofst, 32) ! Length of GRIB message
120  iofst = iofst + 32
121  lensec0 = 16
122  ipos = istart + lensec0
123 
124  ! Currently handles only GRIB Edition 2.
125  if (listsec0(2) .ne. 2) then
126  print *, 'gettemplates: can only decode GRIB edition 2.'
127  ierr = 2
128  return
129  endif
130 
131  ! Loop through the remaining sections keeping track of the length of
132  ! each. Also keep the latest Grid Definition Section info. Unpack
133  ! the requested field number.
134  do
135  ! Check to see if we are at end of GRIB message.
136  ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // cgrib(ipos + 3)
137  if (ctemp .eq. c7777 ) then
138  ipos = ipos+4
139  ! If end of GRIB message not where expected, issue error
140  if (ipos .ne. (istart + lengrib)) then
141  print *, 'gettemplates: "7777" found, but not where expected.'
142  ierr = 4
143  return
144  endif
145  exit
146  endif
147  ! Get length of Section and Section number.
148  iofst = (ipos - 1) * 8
149  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
150  iofst = iofst + 32
151  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Get Section number
152  iofst = iofst + 8
153  !print *, ' lensec = ', lensec, ' secnum = ', isecnum
154 
155  ! If found Section 3, unpack the GDS info using the appropriate
156  ! template. Save in case this is the latest grid before the
157  ! requested field.
158  if (isecnum .eq. 3) then
159  iofst = iofst - 40 ! reset offset to beginning of section
160  call unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, igdslen, ideflist, idefnum, jerr)
161  if (jerr .eq. 0) then
162  have3 = .true.
163  else
164  ierr = 10
165  return
166  endif
167  endif
168 
169  ! If found Section 4, check to see if this field is the one
170  ! requested.
171  if (isecnum .eq. 4) then
172  numfld = numfld + 1
173  if (numfld .eq. ifldnum) then
174  iofst = iofst - 40 ! reset offset to beginning of section
175  call unpack4(cgrib, lcgrib, iofst, ipdsnum, ipdstmpl, ipdslen, coordlist, numcoord, jerr)
176  if (jerr .eq. 0) then
177  have4 = .true.
178  else
179  ierr = 11
180  return
181  endif
182  endif
183  endif
184 
185  ! Check to see if we read pass the end of the GRIB message and
186  ! missed the terminator string '7777'.
187  ipos = ipos + lensec ! Update beginning of section pointer
188  if (ipos .gt. (istart + lengrib)) then
189  print *, 'gettemplates: "7777" not found at end of GRIB message.'
190  ierr = 7
191  return
192  endif
193 
194  if (have3.and.have4) return
195  enddo
196 
197  ! If exited from above loop, the end of the GRIB message was reached
198  ! before the requested field was found.
199  print *, 'gettemplates: GRIB message contained ', numfld, ' different fields.'
200  print *, 'gettemplates: The request was for the ', ifldnum, ' field.'
201  ierr = 6
202 
203 end subroutine gettemplates
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 unpack4(cgrib, lcgrib, iofst, ipdsnum, ipdstmpl, mappdslen, coordlist, numcoord, ierr)
This subroutine unpacks Section 4 (Product Definition Section) starting at octet 6 of that Section.
Definition: getfield.F90:506
subroutine unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, mapgridlen, ideflist, idefnum, ierr)
This subroutine unpacks Section 3 (Grid Definition Section) starting at octet 6 of that Section.
Definition: getfield.F90:376
subroutine gettemplates(cgrib, lcgrib, ifldnum, igds, igdstmpl, igdslen, ideflist, idefnum, ipdsnum, ipdstmpl, ipdslen, coordlist, numcoord, ierr)
This subroutine returns the Grid Definition, and Product Definition for a given data field.