NCEPLIBS-bufr  12.1.0
bitmaps.F90
Go to the documentation of this file.
1 
5 
12 subroutine strbtm ( n, lun )
13 
14  use modv_vars, only: mxbtm, mxbtmse
15 
16  use moda_msgcwd
17  use moda_usrint
18  use moda_tables
19  use moda_bitmaps
20 
21  implicit none
22 
23  integer, intent(in) :: n, lun
24  integer node, nodtam, ii, jj, ibfms, lstjpb
25 
26  logical isbtme
27 
28  node = inv( n, lun )
29 
30  if ( tag(node)(1:5) == 'DPRI ' ) then
31  ! Confirm that this is really an entry within a bitmap. Although it's rare, it is possible for a DPRI element
32  ! to appear in a subset definition outside of a bitmap.
33  isbtme = .false.
34  if ( ntamc > 0 ) then
35  nodtam = lstjpb( node, lun, 'SUB' )
36  do ii = 1, ntamc
37  if ( nodtam == inodtamc(ii) ) then
38  do jj = 1, ntco(ii)
39  if ( ( inodtco(ii,jj) >= inode(lun) ) .and. ( inodtco(ii,jj) <= isc(inode(lun)) ) .and. &
40  ( inodtco(ii,jj) < node ) ) then
41  if ( ctco(ii,jj) == '236000' ) then
42  isbtme = .true.
43  else if ( ( ctco(ii,jj) == '235000' ) .or. ( ctco(ii,jj) == '237255' ) ) then
44  isbtme = .false.
45  end if
46  end if
47  end do
48  end if
49  end do
50  end if
51  if ( .not. isbtme ) then
52  linbtm = .false.
53  return
54  endif
55  if ( .not. linbtm ) then
56  ! This is the start of a new bitmap.
57  if ( nbtm >= mxbtm ) call bort('BUFRLIB: STRBTM - MXBTM OVERFLOW')
58  nbtm = nbtm + 1
59  istbtm(nbtm) = n
60  iszbtm(nbtm) = 0
61  nbtmse(nbtm) = 0
62  linbtm = .true.
63  end if
64  iszbtm(nbtm) = iszbtm(nbtm) + 1
65  if ( ibfms(val(n,lun)) == 0 ) then
66  ! This is a "set" (value=0) entry in the bitmap.
67  if ( nbtmse(nbtm) >= mxbtmse ) call bort('BUFRLIB: STRBTM - MXBTMSE OVERFLOW')
68  nbtmse(nbtm) = nbtmse(nbtm) + 1
70  end if
71  else if ( itp(node) > 1 ) then
72  linbtm = .false.
73  end if
74 
75  return
76 end subroutine strbtm
77 
102 recursive subroutine gettagre ( lunit, tagi, ntagi, tagre, ntagre, iret )
103 
104  use modv_vars, only: im8b
105 
106  use moda_usrint
107  use moda_msgcwd
108  use moda_tables
109 
110  implicit none
111 
112  integer, intent(in) :: lunit, ntagi
113  integer, intent(out) :: iret, ntagre
114  integer my_lunit, my_ntagi, lun, il, im, ni, nre, ltre, ii
115 
116  character*(*), intent(in) :: tagi
117  character*(*), intent(out) :: tagre
118  character*10 tagtmp
119 
120  ! Check for I8 integers.
121 
122  if(im8b) then
123  im8b=.false.
124  call x84(lunit,my_lunit,1)
125  call x84(ntagi,my_ntagi,1)
126  call gettagre(my_lunit,tagi,my_ntagi,tagre,ntagre,iret)
127  call x48(ntagre,ntagre,1)
128  call x48(iret,iret,1)
129  im8b=.true.
130  return
131  endif
132 
133  iret = -1
134 
135  ! Get lun from lunit.
136 
137  call status( lunit, lun, il, im )
138  if ( il == 0 ) return
139  if ( inode(lun) /= inv(1,lun) ) return
140 
141  ! Get tagre and ntagre from the (ntagi)th occurrence of tagi.
142 
143  call fstag( lun, tagi, ntagi, 1, ni, iret )
144  if ( iret /= 0 ) return
145  nre = nrfelm(ni,lun)
146  if ( nre > 0 ) then
147  iret = 0
148  tagre = tag(inv(nre,lun))
149  call strsuc( tagre, tagtmp, ltre )
150  ntagre = 0
151  do ii = 1, nre
152  if ( tag(inv(ii,lun))(1:ltre) == tagre(1:ltre) ) then
153  ntagre = ntagre + 1
154  end if
155  end do
156  end if
157 
158  return
159 end subroutine gettagre
160 
175 integer function igetrfel ( n, lun ) result ( iret )
176 
177  use moda_msgcwd
178  use moda_usrint
179  use moda_tables
180  use moda_bitmaps
181  use moda_nrv203
182 
183  implicit none
184 
185  integer, intent(in) :: n, lun
186  integer node, ii, jj, nn, idxta, idn, ntc, nodflw, nodl236, nodbmap, nodrfe, nodnn, nodtam, idxbtm, iemrk, iect, &
187  lstjpb, imrkopr
188 
189  character*(*), parameter :: bort_str_mrkopr = &
190  'BUFRLIB: IGETRFEL - UNABLE TO FIND PREVIOUS ELEMENT REFERENCED BY MARKER OPERATOR '
191  character*128 bort_str
192  character*6 cflwopr, adn30, fxy
193  character*1 tab
194 
195  iret = 0
196 
197  node = inv( n, lun )
198 
199  if ( itp(node) > 1 ) then
200  if ( node == lstnod ) then
201  lstnodct = lstnodct + 1
202  else
203  lstnod = node
204  lstnodct = 1
205  end if
206  ! Does this subset definition contain any Table C operators with an X value of 21 or greater?
207  idxta = 0
208  if ( ntamc > 0 ) then
209  nodtam = lstjpb( node, lun, 'SUB' )
210  do ii = 1, ntamc
211  if ( nodtam == inodtamc(ii) ) then
212  idxta = ii
213  ntc = ntco(ii)
214  end if
215  end do
216  end if
217  if ( ( idxta > 0 ) .and. ( nbtm > 0 ) ) then
218  ! Check whether this element references a previous element in the same subset via an internal bitmap. To do this,
219  ! we first need to determine the appropriate "follow" operator (if any) corresponding to this element.
220  cflwopr = 'XXXXXX'
221  if ( imrkopr(tag(node)) == 1 ) then
222  cflwopr = tag(node)(1:3) // '000'
223  else
224  call nemtab( lun, tag(node), idn, tab, nn )
225  if ( tab == 'B' ) then
226  fxy = adn30(idn,6)
227  if ( fxy(2:3) == '33' ) cflwopr = '222000'
228  end if
229  end if
230  if ( cflwopr == 'XXXXXX' ) return
231  ! Now, check whether the appropriate "follow" operator was actually present in the subset. If there are multiple
232  ! occurrences, we want the one that most recently precedes the element in question.
233  nodflw = 0
234  do jj = 1, ntc
235  if ( ( ctco(idxta,jj) == cflwopr ) .and. ( inodtco(idxta,jj) >= inode(lun) ) .and. &
236  ( inodtco(idxta,jj) <= isc(inode(lun)) ) .and. ( inodtco(idxta,jj) < node ) ) nodflw = inodtco(idxta,jj)
237  enddo
238  if ( nodflw == 0 ) then
239  if ( imrkopr(tag(node)) == 1 ) then
240  write(bort_str,'("BUFRLB: IGETRFEL - UNABLE TO FIND FOLLOW OPERATOR ",A," IN SUBSET")') cflwopr
241  call bort(bort_str)
242  endif
243  return
244  end if
245  ! We found an appropriate corresponding "follow" operator, so now we need to look for a bitmap corresponding to
246  ! this operator. First, look for a bitmap indicator.
247  nodl236 = 0
248  nodbmap = 0
249  jj = 1
250  do while ( ( jj <= ntc ) .and. ( inodtco(idxta,jj) >= inode(lun) ) .and. &
251  ( inodtco(idxta,jj) <= isc(inode(lun)) ) .and. ( nodbmap == 0 ) )
252  if ( ctco(idxta,jj) == '236000' ) then
253  nodl236 = inodtco(idxta,jj)
254  if ( inodtco(idxta,jj) == nodflw ) nodbmap = nodflw
255  else if ( ( ctco(idxta,jj) == '235000' ) .or. ( ctco(idxta,jj) == '237255' ) ) then
256  nodl236 = 0
257  else if ( ( ctco(idxta,jj) == '237000' ) .and. ( inodtco(idxta,jj) == nodflw ) .and. ( nodl236 /= 0 ) ) then
258  nodbmap = nodl236
259  end if
260  jj = jj + 1
261  end do
262  if ( nodbmap == 0 ) then
263  ! There was no valid bitmap indicator, so we'll just look for a bitmap after the "follow" indicator.
264  nodbmap = nodflw
265  end if
266  ! Find the corresponding bitmap.
267  nn = 1
268  idxbtm = 0
269  do while ( ( idxbtm == 0 ) .and. ( nn <= nval(lun) ) )
270  if ( inv( nn, lun ) > nodbmap ) then
271  ii = 1
272  do while ( ( idxbtm == 0 ) .and. ( ii <= nbtm ) )
273  if ( nn == istbtm(ii) ) then
274  idxbtm = ii
275  else
276  ii = ii + 1
277  end if
278  end do
279  end if
280  nn = nn + 1
281  end do
282  if ( idxbtm == 0 ) then
283  if ( imrkopr(tag(node)) == 1 ) then
284  write(bort_str,'("BUFRLB: IGETRFEL - UNABLE TO FIND BITMAP FOR MARKER OPERATOR ",A)') tag(node)
285  call bort(bort_str)
286  endif
287  return
288  end if
289  ! Use the bitmap to find the previous element in the subset that is referenced by the current element.
290  ! Search backwards from the start of the bitmap, but make sure not to cross a 2-35-000 operator.
291  if ( lstnodct > nbtmse(idxbtm) ) then
292  if ( imrkopr(tag(node)) == 1 ) call bort( bort_str_mrkopr // tag(node) )
293  return
294  end if
295  iemrk = iszbtm(idxbtm) - ibtmse(idxbtm,lstnodct) + 1
296  iect = 0
297  do while ( ( nn >= 1 ) .and. ( iret == 0 ) )
298  nodnn = inv( nn, lun )
299  if ( nodnn <= nodbmap ) then
300  do jj = 1, ntc
301  if ( ( nodnn == inodtco(idxta,jj) ) .and. ( ctco(idxta,jj) == '235000' ) ) then
302  if ( imrkopr(tag(node)) == 1 ) call bort( bort_str_mrkopr // tag(node) )
303  return
304  end if
305  end do
306  if ( itp(nodnn) > 1 ) then
307  iect = iect + 1
308  if ( iect == iemrk ) iret = nn
309  end if
310  end if
311  nn = nn - 1
312  end do
313  if ( iret == 0 ) then
314  if ( imrkopr(tag(node)) == 1 ) call bort( bort_str_mrkopr // tag(node) )
315  return
316  end if
317  if ( imrkopr(tag(node)) == 1 ) then
318  ! This element is a marker operator, so set the scale, reference value and bit width accordingly based on
319  ! those of the previous referenced element.
320  nodrfe = inv( iret, lun )
321  isc(node) = isc(nodrfe)
322  if ( tag(node)(1:3) == '225' ) then
323  ibt(node) = ibt(nodrfe) + 1
324  irf(node) = -1 * (2 ** ibt(nodrfe))
325  else
326  ibt(node) = ibt(nodrfe)
327  irf(node) = irf(nodrfe)
328  if ( nnrv > 0 ) then
329  do ii = 1, nnrv
330  if ( ( nodrfe /= inodnrv(ii) ) .and. ( tag(nodrfe)(1:8) == tagnrv(ii) ) .and. &
331  ( nodrfe >= isnrv(ii) ) .and. ( nodrfe <= ienrv(ii) ) ) then
332  irf(node) = int(nrv(ii))
333  return
334  end if
335  end do
336  end if
337  end if
338  end if
339  end if
340  end if
341 
342  return
343 end function igetrfel
344 
353 integer function imrkopr(nemo) result(iret)
354 
355  implicit none
356 
357  character*(*), intent(in) :: nemo
358 
359  if (len(nemo)<6) then
360  iret = 0
361  else if ( ( nemo(4:6)=='255' ) .and. &
362  ( ( nemo(1:3)=='223' ) .or. ( nemo(1:3)=='224' ) .or. ( nemo(1:3)=='225' ) .or. ( nemo(1:3)=='232' ) ) ) then
363  iret = 1
364  else
365  iret = 0
366  endif
367 
368  return
369 end function imrkopr
integer function igetrfel(n, lun)
Check whether a subset element refers to a previous element within the same subset via an internal bi...
Definition: bitmaps.F90:176
subroutine strbtm(n, lun)
Store internal information in module moda_bitmaps if the input element is part of a bitmap.
Definition: bitmaps.F90:13
recursive subroutine gettagre(lunit, tagi, ntagi, tagre, ntagre, iret)
Check whether a specified Table B mnemonic references another Table B mnemonic within the same data s...
Definition: bitmaps.F90:103
integer function imrkopr(nemo)
Check whether a specified mnemonic is a Table C marker operator.
Definition: bitmaps.F90:354
subroutine bort(str)
Log an error message, then abort the application program.
Definition: borts.F90:15
subroutine nemtab(lun, nemo, idn, tab, iret)
Get information about a descriptor, based on a mnemonic.
Definition: fxy.F90:432
character *(*) function adn30(idn, ldn)
Convert an FXY value from its WMO bit-wise representation to a character string of length 5 or 6.
Definition: fxy.F90:18
subroutine strsuc(str1, str2, lens)
Remove leading and trailing blanks from a character string.
Definition: misc.F90:220
Declare arrays and variables used to store bitmaps internally within a data subset definition.
integer, dimension(:), allocatable iszbtm
Size of bitmap (total number of entries, whether "set" (set to a value of 0) or not).
integer lstnod
Most recent jump/link table entry that was processed by function igetrfel() and whose corresponding v...
integer, dimension(:,:), allocatable inodtco
Entries within jump/link table which contain Table C operators.
integer ntamc
Number of Table A mnemonics in jump/link table (up to a maximum of mxtamc) which contain at least one...
integer, dimension(:), allocatable istbtm
Ordinal position in data subset definition corresponding to the first entry of the bitmap.
integer nbtm
Number of stored bitmaps for the current data subset (up to a maximum of mxbtm).
integer, dimension(:), allocatable inodtamc
Entries within jump/link table which contain Table A mnemonics.
integer, dimension(:,:), allocatable ibtmse
Ordinal positions in bitmap of bits that were "set" (set to a value of 0); these ordinal positions ca...
integer lstnodct
Current count of consecutive occurrences of lstnod.
logical linbtm
true if a bitmap is in the process of being read for the current data subset; false otherwise.
integer, dimension(:), allocatable nbtmse
Number of "set" entries (set to a value of 0) in the bitmap.
character *6, dimension(:,:), allocatable ctco
Table C operators corresponding to inodtco.
integer, dimension(:), allocatable ntco
Number of Table C operators (with an XX value of 21 or greater) within the data subset definition of ...
Declare arrays used to store information about the current BUFR message that is in the process of bei...
integer, dimension(:), allocatable inode
Table A mnemonic for type of BUFR message.
Declare arrays and variables for use with any 2-03-YYY (change reference value) operators present wit...
integer, dimension(:), allocatable ienrv
End of entry range in jump/link table, within which the corresponding new reference value in nrv will...
character *8, dimension(:), allocatable tagnrv
Table B mnemonic to which the corresponding new reference value in nrv applies.
integer, dimension(:), allocatable isnrv
Start of entry range in jump/link table, within which the corresponding new reference value in nrv wi...
integer nnrv
Number of entries in the jump/link table which contain new reference values (up to a maximum of mxnrv...
integer *8, dimension(:), allocatable nrv
New reference values corresponding to inodnrv.
integer, dimension(:), allocatable inodnrv
Entries within jump/link table which contain new reference values.
Declare arrays and variables used to store the internal jump/link table.
integer, dimension(:), allocatable irf
Reference values corresponding to tag and typ:
integer, dimension(:), allocatable isc
Scale factors corresponding to tag and typ:
integer, dimension(:), allocatable ibt
Bit widths corresponding to tag and typ:
character *10, dimension(:), allocatable tag
Mnemonics in the jump/link table.
integer, dimension(:), allocatable itp
Integer type values corresponding to typ:
Declare arrays used to store data values and associated metadata for the current BUFR data subset in ...
integer, dimension(:), allocatable nval
Number of data values in BUFR data subset.
real *8, dimension(:,:), allocatable, target val
Data values.
integer, dimension(:,:), allocatable, target inv
Inventory pointer which links each data value to its corresponding node in the internal jump/link tab...
integer, dimension(:,:), allocatable nrfelm
Referenced data value, for data values which refer to a previous data value in the BUFR data subset v...
recursive subroutine status(lunit, lun, il, im)
Check whether a specified Fortran logical unit number is currently connected to the NCEPLIBS-bufr sof...
subroutine x48(iin4, iout8, nval)
Encode one or more 4-byte integer values as 8-byte integer values.
Definition: x4884.F90:18
subroutine x84(iin8, iout4, nval)
Encode one or more 8-byte integer values as 4-byte integer values.
Definition: x4884.F90:65