NCEPLIBS-bufr  12.3.0
bitmaps.F90
Go to the documentation of this file.
1 
5 
19 subroutine strbtm ( n, lun, ival )
20 
21  use modv_vars, only: mxbtm, mxbtmse
22 
23  use moda_msgcwd
24  use moda_usrint
25  use moda_tables
26  use moda_bitmaps
27 
28  implicit none
29 
30  integer, intent(in) :: n, lun, ival
31  integer node, nodtam, ii, jj, lstjpb
32 
33  logical isbtme
34 
35  node = inv( n, lun )
36 
37  if ( tag(node)(1:5) == 'DPRI ' ) then
38  ! Confirm that this is really an entry within a bitmap. Although it's rare, it is possible for a DPRI element
39  ! to appear in a subset definition outside of a bitmap.
40  isbtme = .false.
41  if ( ntamc > 0 ) then
42  nodtam = lstjpb( node, lun, 'SUB' )
43  do ii = 1, ntamc
44  if ( nodtam == inodtamc(ii) ) then
45  do jj = 1, ntco(ii)
46  if ( ( inodtco(ii,jj) >= inode(lun) ) .and. ( inodtco(ii,jj) <= isc(inode(lun)) ) .and. &
47  ( inodtco(ii,jj) < node ) ) then
48  if ( ctco(ii,jj) == '236000' ) then
49  isbtme = .true.
50  else if ( ( ctco(ii,jj) == '235000' ) .or. ( ctco(ii,jj) == '237255' ) ) then
51  isbtme = .false.
52  end if
53  end if
54  end do
55  end if
56  end do
57  end if
58  if ( .not. isbtme ) then
59  linbtm = .false.
60  return
61  endif
62  if ( .not. linbtm ) then
63  ! This is the start of a new bitmap.
64  if ( nbtm >= mxbtm ) call bort('BUFRLIB: STRBTM - MXBTM OVERFLOW')
65  nbtm = nbtm + 1
66  istbtm(nbtm) = n
67  iszbtm(nbtm) = 0
68  nbtmse(nbtm) = 0
69  linbtm = .true.
70  end if
71  iszbtm(nbtm) = iszbtm(nbtm) + 1
72  if ( ival == 0 ) then
73  ! This is a "set" (value=0) entry in the bitmap.
74  if ( nbtmse(nbtm) >= mxbtmse ) call bort('BUFRLIB: STRBTM - MXBTMSE OVERFLOW')
75  nbtmse(nbtm) = nbtmse(nbtm) + 1
77  end if
78  else if ( itp(node) > 1 ) then
79  linbtm = .false.
80  end if
81 
82  return
83 end subroutine strbtm
84 
109 recursive subroutine gettagre ( lunit, tagi, ntagi, tagre, ntagre, iret )
110 
111  use bufrlib
112 
113  use modv_vars, only: im8b
114 
115  use moda_usrint
116  use moda_msgcwd
117  use moda_tables
118 
119  implicit none
120 
121  integer, intent(in) :: lunit, ntagi
122  integer, intent(out) :: iret, ntagre
123  integer my_lunit, my_ntagi, lun, il, im, ni, nre, ltre, ii, lci, ntrchr, ltr, bort_target_set
124 
125  character*(*), intent(in) :: tagi
126  character*(*), intent(out) :: tagre
127  character*10 tagtmp, ctagi, ctagre
128 
129  ! Check for I8 integers.
130 
131  if(im8b) then
132  im8b=.false.
133  call x84(lunit,my_lunit,1)
134  call x84(ntagi,my_ntagi,1)
135  call gettagre(my_lunit,tagi,my_ntagi,tagre,ntagre,iret)
136  call x48(ntagre,ntagre,1)
137  call x48(iret,iret,1)
138  im8b=.true.
139  return
140  endif
141 
142  tagre = ' '
143 
144  ! If we're catching bort errors, set a target return location if one doesn't already exist.
145 
146  if ( bort_target_set() == 1 ) then
147  call strsuc( tagi, ctagi, lci )
148  call catch_bort_gettagre_c( lunit, ctagi, lci, ntagi, ctagre, len(ctagre), ntagre, ntrchr, iret )
149  ltr = min( len(tagre), ntrchr )
150  tagre(1:ltr) = ctagre(1:ltr)
151  call bort_target_unset
152  return
153  endif
154 
155  iret = -1
156 
157  ! Get lun from lunit.
158 
159  call status( lunit, lun, il, im )
160  if ( il == 0 ) return
161  if ( inode(lun) /= inv(1,lun) ) return
162 
163  ! Get tagre and ntagre from the (ntagi)th occurrence of tagi.
164 
165  call fstag( lun, tagi, ntagi, 1, ni, iret )
166  if ( iret /= 0 ) return
167  nre = nrfelm(ni,lun)
168  if ( nre > 0 ) then
169  iret = 0
170  tagre = tag(inv(nre,lun))
171  call strsuc( tagre, tagtmp, ltre )
172  ntagre = 0
173  do ii = 1, nre
174  if ( tag(inv(ii,lun))(1:ltre) == tagre(1:ltre) ) then
175  ntagre = ntagre + 1
176  end if
177  end do
178  end if
179 
180  return
181 end subroutine gettagre
182 
197 integer function igetrfel ( n, lun ) result ( iret )
198 
199  use moda_msgcwd
200  use moda_usrint
201  use moda_tables
202  use moda_bitmaps
203  use moda_nrv203
204 
205  implicit none
206 
207  integer, intent(in) :: n, lun
208  integer node, ii, jj, nn, idxta, idn, ntc, nodflw, nodl236, nodbmap, nodrfe, nodnn, nodtam, idxbtm, iemrk, iect, &
209  lstjpb, imrkopr
210 
211  character*(*), parameter :: bort_str_mrkopr = &
212  'BUFRLIB: IGETRFEL - UNABLE TO FIND PREVIOUS ELEMENT REFERENCED BY MARKER OPERATOR '
213  character*128 bort_str
214  character*6 cflwopr, adn30, fxy
215  character*1 tab
216 
217  iret = 0
218 
219  node = inv( n, lun )
220 
221  if ( itp(node) > 1 ) then
222  if ( node == lstnod ) then
223  lstnodct = lstnodct + 1
224  else
225  lstnod = node
226  lstnodct = 1
227  end if
228  ! Does this subset definition contain any Table C operators with an X value of 21 or greater?
229  idxta = 0
230  if ( ntamc > 0 ) then
231  nodtam = lstjpb( node, lun, 'SUB' )
232  do ii = 1, ntamc
233  if ( nodtam == inodtamc(ii) ) then
234  idxta = ii
235  ntc = ntco(ii)
236  end if
237  end do
238  end if
239  if ( ( idxta > 0 ) .and. ( nbtm > 0 ) ) then
240  ! Check whether this element references a previous element in the same subset via an internal bitmap. To do this,
241  ! we first need to determine the appropriate "follow" operator (if any) corresponding to this element.
242  cflwopr = 'XXXXXX'
243  if ( imrkopr(tag(node)) == 1 ) then
244  cflwopr = tag(node)(1:3) // '000'
245  else
246  call nemtab( lun, tag(node), idn, tab, nn )
247  if ( tab == 'B' ) then
248  fxy = adn30(idn,6)
249  if ( fxy(2:3) == '33' ) cflwopr = '222000'
250  end if
251  end if
252  if ( cflwopr == 'XXXXXX' ) return
253  ! Now, check whether the appropriate "follow" operator was actually present in the subset. If there are multiple
254  ! occurrences, we want the one that most recently precedes the element in question.
255  nodflw = 0
256  do jj = 1, ntc
257  if ( ( ctco(idxta,jj) == cflwopr ) .and. ( inodtco(idxta,jj) >= inode(lun) ) .and. &
258  ( inodtco(idxta,jj) <= isc(inode(lun)) ) .and. ( inodtco(idxta,jj) < node ) ) nodflw = inodtco(idxta,jj)
259  enddo
260  if ( nodflw == 0 ) then
261  if ( imrkopr(tag(node)) == 1 ) then
262  write(bort_str,'("BUFRLB: IGETRFEL - UNABLE TO FIND FOLLOW OPERATOR ",A," IN SUBSET")') cflwopr
263  call bort(bort_str)
264  endif
265  return
266  end if
267  ! We found an appropriate corresponding "follow" operator, so now we need to look for a bitmap corresponding to
268  ! this operator. First, look for a bitmap indicator.
269  nodl236 = 0
270  nodbmap = 0
271  jj = 1
272  do while ( ( jj <= ntc ) .and. ( inodtco(idxta,jj) >= inode(lun) ) .and. &
273  ( inodtco(idxta,jj) <= isc(inode(lun)) ) .and. ( nodbmap == 0 ) )
274  if ( ctco(idxta,jj) == '236000' ) then
275  nodl236 = inodtco(idxta,jj)
276  if ( inodtco(idxta,jj) == nodflw ) nodbmap = nodflw
277  else if ( ( ctco(idxta,jj) == '235000' ) .or. ( ctco(idxta,jj) == '237255' ) ) then
278  nodl236 = 0
279  else if ( ( ctco(idxta,jj) == '237000' ) .and. ( inodtco(idxta,jj) == nodflw ) .and. ( nodl236 /= 0 ) ) then
280  nodbmap = nodl236
281  end if
282  jj = jj + 1
283  end do
284  if ( nodbmap == 0 ) then
285  ! There was no valid bitmap indicator, so we'll just look for a bitmap after the "follow" indicator.
286  nodbmap = nodflw
287  end if
288  ! Find the corresponding bitmap.
289  nn = 1
290  idxbtm = 0
291  do while ( ( idxbtm == 0 ) .and. ( nn <= nval(lun) ) )
292  if ( inv( nn, lun ) > nodbmap ) then
293  ii = 1
294  do while ( ( idxbtm == 0 ) .and. ( ii <= nbtm ) )
295  if ( nn == istbtm(ii) ) then
296  idxbtm = ii
297  else
298  ii = ii + 1
299  end if
300  end do
301  end if
302  nn = nn + 1
303  end do
304  if ( idxbtm == 0 ) then
305  if ( imrkopr(tag(node)) == 1 ) then
306  write(bort_str,'("BUFRLB: IGETRFEL - UNABLE TO FIND BITMAP FOR MARKER OPERATOR ",A)') tag(node)
307  call bort(bort_str)
308  endif
309  return
310  end if
311  ! Use the bitmap to find the previous element in the subset that is referenced by the current element.
312  ! Search backwards from the start of the bitmap, but make sure not to cross a 2-35-000 operator.
313  if ( lstnodct > nbtmse(idxbtm) ) then
314  if ( imrkopr(tag(node)) == 1 ) call bort( bort_str_mrkopr // tag(node) )
315  return
316  end if
317  iemrk = iszbtm(idxbtm) - ibtmse(idxbtm,lstnodct) + 1
318  iect = 0
319  do while ( ( nn >= 1 ) .and. ( iret == 0 ) )
320  nodnn = inv( nn, lun )
321  if ( nodnn <= nodbmap ) then
322  do jj = 1, ntc
323  if ( ( nodnn == inodtco(idxta,jj) ) .and. ( ctco(idxta,jj) == '235000' ) ) then
324  if ( imrkopr(tag(node)) == 1 ) call bort( bort_str_mrkopr // tag(node) )
325  return
326  end if
327  end do
328  if ( itp(nodnn) > 1 ) then
329  iect = iect + 1
330  if ( iect == iemrk ) iret = nn
331  end if
332  end if
333  nn = nn - 1
334  end do
335  if ( iret == 0 ) then
336  if ( imrkopr(tag(node)) == 1 ) call bort( bort_str_mrkopr // tag(node) )
337  return
338  end if
339  if ( imrkopr(tag(node)) == 1 ) then
340  ! This element is a marker operator, so set the scale, reference value and bit width accordingly based on
341  ! those of the previous referenced element.
342  nodrfe = inv( iret, lun )
343  isc(node) = isc(nodrfe)
344  if ( tag(node)(1:3) == '225' ) then
345  ibt(node) = ibt(nodrfe) + 1
346  irf(node) = -1 * (2 ** ibt(nodrfe))
347  else
348  ibt(node) = ibt(nodrfe)
349  irf(node) = irf(nodrfe)
350  if ( nnrv > 0 ) then
351  do ii = 1, nnrv
352  if ( ( nodrfe /= inodnrv(ii) ) .and. ( tag(nodrfe)(1:8) == tagnrv(ii) ) .and. &
353  ( nodrfe >= isnrv(ii) ) .and. ( nodrfe <= ienrv(ii) ) ) then
354  irf(node) = int(nrv(ii))
355  return
356  end if
357  end do
358  end if
359  end if
360  end if
361  end if
362  end if
363 
364  return
365 end function igetrfel
366 
375 integer function imrkopr(nemo) result(iret)
376 
377  implicit none
378 
379  character*(*), intent(in) :: nemo
380 
381  if (len(nemo)<6) then
382  iret = 0
383  else if ( ( nemo(4:6)=='255' ) .and. &
384  ( ( nemo(1:3)=='223' ) .or. ( nemo(1:3)=='224' ) .or. ( nemo(1:3)=='225' ) .or. ( nemo(1:3)=='232' ) ) ) then
385  iret = 1
386  else
387  iret = 0
388  endif
389 
390  return
391 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:198
subroutine strbtm(n, lun, ival)
Store internal information in module moda_bitmaps if the input element is part of a bitmap.
Definition: bitmaps.F90:20
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:110
integer function imrkopr(nemo)
Check whether a specified mnemonic is a Table C marker operator.
Definition: bitmaps.F90:376
recursive subroutine bort(str)
Log an error message, then either return to or abort the application program.
Definition: borts.F90:15
subroutine bort_target_unset
Clear any existing bort target.
Definition: borts.F90:180
integer function bort_target_set()
Sets a new bort target, if bort catching is enabled and such a target doesn't already exist.
Definition: borts.F90:160
subroutine nemtab(lun, nemo, idn, tab, iret)
Get information about a descriptor, based on a mnemonic.
Definition: fxy.F90:434
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:199
Wrap C NCEPLIBS-bufr functions so they can be called from within the Fortran part of the library.
Definition: bufrlib.F90:11
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
Set to .true.
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.
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