NCEPLIBS-bufr  12.2.0
All Data Structures Namespaces Files Functions Variables Macros Pages
readwritesb.F90
Go to the documentation of this file.
1 
5 
31 recursive subroutine readsb(lunit,iret)
32 
33  use modv_vars, only: im8b
34 
35  use moda_msgcwd
36  use moda_unptyp
37  use moda_bitbuf
38  use moda_bitmaps
39  use moda_stcode
40 
41  implicit none
42 
43  integer, intent(in) :: lunit
44  integer, intent(out) :: iret
45  integer my_lunit, lun, il, im, ier, nbyt
46 
47  ! Check for I8 integers
48 
49  if(im8b) then
50  im8b=.false.
51 
52  call x84(lunit,my_lunit,1)
53  call readsb(my_lunit,iret)
54  call x48(iret,iret,1)
55 
56  im8b=.true.
57  return
58  endif
59 
60  iret = -1
61 
62  ! Check the file status
63 
64  call status(lunit,lun,il,im)
65  if(il==0) call bort('BUFRLIB: READSB - INPUT BUFR FILE IS CLOSED, IT MUST BE OPEN FOR INPUT')
66  if(il>0) call bort('BUFRLIB: READSB - INPUT BUFR FILE IS OPEN FOR OUTPUT, IT MUST BE OPEN FOR INPUT')
67  if(im==0) return
68 
69  ! See if there is another subset in the message
70 
71  if(nsub(lun)==msub(lun)) return
72  nsub(lun) = nsub(lun) + 1
73 
74  ! Read the next subset and reset the pointers
75 
76  nbtm = 0
77  lstnod = 0
78  lstnodct = 0
79  iscodes(lun) = 0
80  linbtm = .false.
81 
82  if(msgunp(lun)==0) then
83  ibit = mbyt(lun)*8
84  call upb(nbyt,16,mbay(1,lun),ibit)
85  call rdtree(lun,ier)
86  if(ier/=0) return
87  mbyt(lun) = mbyt(lun) + nbyt
88  elseif(msgunp(lun)==1) then
89  ! message with "standard" Section 3
90  ibit = mbyt(lun)
91  call rdtree(lun,ier)
92  if(ier/=0) return
93  mbyt(lun) = ibit
94  else
95  ! compressed message
96  call rdcmps(lun)
97  if (iscodes(lun) /= 0) return
98  endif
99 
100  iret = 0
101 
102  return
103 end subroutine readsb
104 
116 recursive integer function ireadsb(lunit) result(iret)
117 
118  use modv_vars, only: im8b
119 
120  implicit none
121 
122  integer, intent(in) :: lunit
123  integer my_lunit
124 
125  ! Check for I8 integers.
126 
127  if(im8b) then
128  im8b=.false.
129 
130  call x84(lunit,my_lunit,1)
131  iret=ireadsb(my_lunit)
132 
133  im8b=.true.
134  return
135  endif
136 
137  call readsb(lunit,iret)
138 
139  return
140 end function ireadsb
141 
167 recursive subroutine readns(lunit,subset,jdate,iret)
168 
169  use modv_vars, only: im8b
170 
171  use moda_msgcwd
172  use moda_tables
173 
174  implicit none
175 
176  integer, intent(in) :: lunit
177  integer, intent(out) :: jdate, iret
178  integer my_lunit, lun, il, im
179 
180  character*8, intent(out) :: subset
181 
182  ! Check for I8 integers
183 
184  if(im8b) then
185  im8b=.false.
186 
187  call x84(lunit,my_lunit,1)
188  call readns(my_lunit,subset,jdate,iret)
189  call x48(jdate,jdate,1)
190  call x48(iret,iret,1)
191 
192  im8b=.true.
193  return
194  endif
195 
196  ! Refresh the subset and jdate parameters
197 
198  call status(lunit,lun,il,im)
199  if(il==0) call bort('BUFRLIB: READNS - INPUT BUFR FILE IS CLOSED, IT MUST BE OPEN FOR INPUT')
200  if(il>0) call bort('BUFRLIB: READNS - INPUT BUFR FILE IS OPEN FOR OUTPUT, IT MUST BE OPEN FOR INPUT')
201  if(inode(lun)==0) then
202  subset = ' '
203  else
204  subset = tag(inode(lun))(1:8)
205  endif
206  jdate = idate(lun)
207 
208  ! Read the next subset in the BUFR file
209 
210  do while (.true.)
211  call readsb(lunit,iret)
212  if (iret==0) exit
213  call readmg(lunit,subset,jdate,iret)
214  if (iret/=0) exit
215  enddo
216 
217  return
218 end subroutine readns
219 
238 recursive integer function ireadns(lunit,subset,idate) result(iret)
239 
240  use modv_vars, only: im8b
241 
242  implicit none
243 
244  integer, intent(in) :: lunit
245  integer, intent(out) :: idate
246  integer my_lunit
247 
248  character*8, intent(out) :: subset
249 
250  ! Check for I8 integers.
251 
252  if(im8b) then
253  im8b=.false.
254 
255  call x84(lunit,my_lunit,1)
256  iret=ireadns(my_lunit,subset,idate)
257  call x48(idate,idate,1)
258 
259  im8b=.true.
260  return
261  endif
262 
263  call readns(lunit,subset,idate,iret)
264 
265  return
266 end function ireadns
267 
302 recursive subroutine writsb(lunit)
303 
304  use modv_vars, only: im8b
305 
306  use moda_msgcmp
307 
308  implicit none
309 
310  integer, intent(in) :: lunit
311  integer my_lunit, lun, il, im
312 
313  ! Check for I8 integers
314 
315  if(im8b) then
316  im8b=.false.
317 
318  call x84 ( lunit, my_lunit, 1 )
319  call writsb ( my_lunit )
320 
321  im8b=.true.
322  return
323  endif
324 
325  ! Check the file status
326 
327  call status(lunit,lun,il,im)
328  if(il==0) call bort('BUFRLIB: WRITSB - OUTPUT BUFR FILE IS CLOSED, IT MUST BE OPEN FOR OUTPUT')
329  if(il<0) call bort('BUFRLIB: WRITSB - OUTPUT BUFR FILE IS OPEN FOR INPUT, IT MUST BE OPEN FOR OUTPUT')
330  if(im==0) call bort('BUFRLIB: WRITSB - A MESSAGE MUST BE OPEN IN OUTPUT BUFR FILE, NONE ARE')
331 
332  ! Pack up the subset and put it into the message
333 
334  call wrtree(lun)
335  if( ccmf=='Y' ) then
336  call wrcmps(lunit)
337  else
338  call msgupd(lunit,lun)
339  endif
340 
341  return
342 end subroutine writsb
343 
422 recursive subroutine writsa(lunxx,lmsgt,msgt,msgl)
423 
424  use modv_vars, only: im8b
425 
426  use moda_bufrmg
427  use moda_msgcmp
428 
429  implicit none
430 
431  integer, intent(in) :: lunxx, lmsgt
432  integer, intent(out) :: msgt(*), msgl
433  integer my_lunxx, my_lmsgt, lunit, lun, il, im, n
434 
435  ! Check for I8 integers
436 
437  if(im8b) then
438  im8b=.false.
439 
440  call x84 ( lunxx, my_lunxx, 1 )
441  call x84 ( lmsgt, my_lmsgt, 1 )
442  call writsa ( my_lunxx, my_lmsgt*2, msgt, msgl )
443  msgl = msgl/2
444  call x48 ( msgl, msgl, 1 )
445 
446  im8b=.true.
447  return
448  endif
449 
450  lunit = abs(lunxx)
451 
452  ! Check the file status
453 
454  call status(lunit,lun,il,im)
455  if(il==0) call bort('BUFRLIB: WRITSA - OUTPUT BUFR FILE IS CLOSED, IT MUST BE OPEN FOR OUTPUT')
456  if(il<0) call bort('BUFRLIB: WRITSA - OUTPUT BUFR FILE IS OPEN FOR INPUT, IT MUST BE OPEN FOR OUTPUT')
457  if(im==0) call bort('BUFRLIB: WRITSA - A MESSAGE MUST BE OPEN IN OUTPUT BUFR FILE, NONE ARE')
458 
459  ! If lunxx < 0, force memory msg to be written (w/o any current subset)
460 
461  if(lunxx<0) call closmg(lunit)
462 
463  ! Is there a completed BUFR message to be returned?
464 
465  if(msglen(lun)>0) then
466  if(msglen(lun)>lmsgt) call bort('BUFRLIB: WRITSA - OVERFLOW OF OUTPUT BUFR MESSAGE ARRAY; TRY A LARGER '// &
467  'DIMENSION FOR THIS ARRAY')
468  msgl = msglen(lun)
469  do n=1,msgl
470  msgt(n) = msgtxt(n,lun)
471  enddo
472  msglen(lun) = 0
473  else
474  msgl = 0
475  endif
476 
477  if(lunxx<0) return
478 
479  ! Pack up the subset and put it into the message
480 
481  call wrtree(lun)
482  if( ccmf=='Y' ) then
483  call wrcmps(lunit)
484  else
485  call msgupd(lunit,lun)
486  endif
487 
488  ! If the just-completed call to wrcmps() or msgupd() for this subset caused a message to be flushed to abs(lunxx), then
489  ! attempt to retrieve and return that message now. Otherwise, we run the risk that the next call to openmb() or openmg()
490  ! might cause another message to be flushed, and thus overwrite the current message within array msgtxt before we
491  ! had the chance to retrieve it during the next call to writsa().
492 
493  ! Also note that, in rare instances (e.g. if the byte count of the most recent subset is > 65530), we could end up with
494  ! two BUFR messages available to be returned from this one call to writsa(). If sufficient space is available in the
495  ! msgt array, then go ahead and return both messages now.
496 
497  if( (msglen(lun)>0) .and. (msgl+msglen(lun)<=lmsgt) ) then
498  do n = 1,msglen(lun)
499  msgt(msgl+n) = msgtxt(n,lun)
500  enddo
501  msgl = msgl+msglen(lun)
502  msglen(lun) = 0
503  endif
504 
505  return
506 end subroutine writsa
507 
535 recursive subroutine rdmgsb(lunit,imsg,isub)
536 
537  use modv_vars, only: im8b
538 
539  use moda_msgcwd
540  use moda_bitbuf
541 
542  implicit none
543 
544  integer, intent(in) :: lunit, imsg, isub
545  integer my_lunit, my_imsg, my_isub, lun, il, im, i, jdate, iret
546 
547  character*128 bort_str
548  character*8 subset
549 
550  ! Check for I8 integers
551 
552  if(im8b) then
553  im8b=.false.
554 
555  call x84(lunit,my_lunit,1)
556  call x84(imsg,my_imsg,1)
557  call x84(isub,my_isub,1)
558  call rdmgsb(my_lunit,my_imsg,my_isub)
559 
560  im8b=.true.
561  return
562  endif
563 
564  ! Open the file and skip to message #imsg
565 
566  call openbf(lunit,'IN',lunit)
567  call status(lunit,lun,il,im)
568 
569  ! Note that we need to use subroutine readmg() to actually read in all of the messages (including the
570  ! first (imsg-1) messages!), just in case there are any embedded dictionary messages in the file.
571 
572  do i=1,imsg
573  call readmg(lunit,subset,jdate,iret)
574  if(iret<0) then
575  write(bort_str,'("BUFRLIB: RDMGSB - HIT END OF FILE BEFORE READING REQUESTED MESSAGE NO.",I5," IN '//&
576  'BUFR FILE CONNECTED TO UNIT",I4)') imsg,lunit
577  call bort(bort_str)
578  endif
579  enddo
580 
581  ! Position at subset #isub
582 
583  do i=1,isub
584  call readsb(lunit,iret)
585  if(iret<0) then
586  write(bort_str,'("BUFRLIB: RDMGSB - ALL SUBSETS READ BEFORE READING REQ. SUBSET NO.",I3," IN '// &
587  'REQ. MSG NO.",I5," IN BUFR FILE CONNECTED TO UNIT",I4)') isub,imsg,lunit
588  call bort(bort_str)
589  endif
590  enddo
591 
592  return
593 end subroutine rdmgsb
594 
612 subroutine msgupd(lunit,lun)
613 
614  use modv_vars, only: iprt, nby0, nby1, nby2, nby3
615 
616  use moda_msgcwd
617  use moda_bitbuf
618  use moda_h4wlc
619 
620  implicit none
621 
622  integer, intent(in) :: lunit, lun
623  integer ibyt, lbyt, lbit, nbyt, ii, iupb
624 
625  logical msgfull
626 
627  character*128 errstr
628 
629  ! Pad the subset buffer
630 
631  call pad(ibay,ibit,ibyt,8)
632 
633  ! Check whether the new subset should be written into the currently open message
634 
635  if(msgfull(mbyt(lun),ibyt,maxbyt) .or. ((ibyt>65530).and.(nsub(lun)>0))) then
636  ! No it should not, either because it doesn't fit
637  ! OR
638  ! It has byte count > 65530 (sufficiently close to the upper limit for the 16 bit byte counter placed at the beginning
639  ! of each subset), and the current message has at least one subset in it
640  !
641  ! In either of these cases, we need to write out the current message and then create a new one to hold the current subset
642  call msgwrt(lunit,mbay(1,lun),mbyt(lun))
643  call msgini(lun)
644  endif
645 
646  if(msgfull(mbyt(lun),ibyt,maxbyt)) then
647  ! This is an overlarge subset that won't fit in any message given the current value of maxbyt, so discard the subset
648  ! and exit gracefully.
649  if(iprt>=0) then
650  call errwrt('+++++++++++++++++++++WARNING+++++++++++++++++++++++')
651  write ( unit=errstr, fmt='(A,A,I7,A)') 'BUFRLIB: MSGUPD - SUBSET LONGER THAN ANY POSSIBLE MESSAGE ', &
652  '{MAXIMUM MESSAGE LENGTH = ', maxbyt, '}'
653  call errwrt(errstr)
654  call errwrt('>>>>>>>OVERLARGE SUBSET DISCARDED FROM FILE<<<<<<<<')
655  call errwrt('+++++++++++++++++++++WARNING+++++++++++++++++++++++')
656  call errwrt(' ')
657  endif
658  call usrtpl(lun,1,1)
659  return
660  endif
661 
662  ! Set a byte count and transfer the subset buffer into the message
663 
664  lbit = 0
665  call pkb(ibyt,16,ibay,lbit)
666 
667  ! Note that we want to append the data for this subset to the end of Section 4, but the value in mbyt(lun) already includes
668  ! the length of Section 5 (i.e. 4 bytes). Therefore, we need to begin writing at the point 3 bytes prior to the byte
669  ! currently pointed to by mbyt(lun).
670 
671  call mvb(ibay,1,mbay(1,lun),mbyt(lun)-3,ibyt)
672 
673  ! Update the subset and byte counters
674 
675  mbyt(lun) = mbyt(lun) + ibyt
676  nsub(lun) = nsub(lun) + 1
677 
678  lbit = (nby0+nby1+nby2+4)*8
679  call pkb(nsub(lun),16,mbay(1,lun),lbit)
680 
681  lbyt = nby0+nby1+nby2+nby3
682  nbyt = iupb(mbay(1,lun),lbyt+1,24)
683  lbit = lbyt*8
684  call pkb(nbyt+ibyt,24,mbay(1,lun),lbit)
685 
686  ! If any long character strings are being held internally for storage into this subset, store them now
687 
688  if(nh4wlc>0) then
689  do ii = 1, nh4wlc
690  call writlc(luh4wlc(ii),chh4wlc(ii),sth4wlc(ii))
691  enddo
692  nh4wlc = 0
693  endif
694 
695  ! If the subset byte count is > 65530, then give it its own one-subset message (cannot have any other subsets in this
696  ! message because their beginning would be beyond the upper limit of 65535 in the 16-bit byte counter, meaning they
697  ! could not be located!)
698 
699  if(ibyt>65530) then
700  if(iprt>=1) then
701  call errwrt('+++++++++++++++++++++WARNING+++++++++++++++++++++++')
702  write ( unit=errstr, fmt='(A,I7,A,A)') 'BUFRLIB: MSGUPD - SUBSET HAS BYTE COUNT = ',ibyt,' > UPPER LIMIT OF 65535'
703  call errwrt(errstr)
704  call errwrt('>>>>>>>WILL BE WRITTEN INTO ITS OWN MESSAGE<<<<<<<<')
705  call errwrt('+++++++++++++++++++++WARNING+++++++++++++++++++++++')
706  call errwrt(' ')
707  endif
708  call msgwrt(lunit,mbay(1,lun),mbyt(lun))
709  call msgini(lun)
710  endif
711 
712  ! Reset the user arrays
713 
714  call usrtpl(lun,1,1)
715 
716  return
717 end subroutine msgupd
718 
748 subroutine pad(ibay,ibit,ibyt,ipadb)
749 
750  implicit none
751 
752  integer, intent(inout) :: ibay(*), ibit
753  integer, intent(in) :: ipadb
754  integer, intent(out) :: ibyt
755  integer ipad
756 
757  character*128 bort_str
758 
759  ! Pad the subset to an ipadb bit boundary
760 
761  ipad = ipadb - mod(ibit+8,ipadb)
762  ! First pack the # of bits being padded (this is a delayed replication factor)
763  call pkb(ipad,8,ibay,ibit)
764  ! Now pad with zeroes to the byte boundary
765  call pkb(0,ipad,ibay,ibit)
766  ibyt = ibit/8
767 
768  if(mod(ibit,8)/=0) then
769  write(bort_str,'("BUFRLIB: PAD - THE NUMBER OF BITS IN A PACKED'// &
770  ' SUBSET AFTER PADDING (",I8,") IS NOT A MULTIPLE OF 8")') ibit
771  call bort(bort_str)
772  endif
773 
774  return
775 end subroutine pad
776 
802 recursive integer function lcmgdf(lunit,subset) result(iret)
803 
804  use modv_vars, only: im8b
805 
806  use moda_tables
807 
808  implicit none
809 
810  integer, intent(in) :: lunit
811  integer my_lunit, lun, il, im, mtyp, msbt, inod, nte, i
812 
813  character*8, intent(in) :: subset
814 
815  ! Check for I8 integers.
816 
817  if(im8b) then
818  im8b=.false.
819 
820  call x84(lunit,my_lunit,1)
821  iret=lcmgdf(my_lunit,subset)
822 
823  im8b=.true.
824  return
825  endif
826 
827  iret = 0
828 
829  ! Get lun from lunit.
830 
831  call status(lunit,lun,il,im)
832  if (il==0) call bort('BUFRLIB: LCMGDF - INPUT BUFR FILE IS CLOSED, IT MUST BE OPEN')
833 
834  ! Confirm that subset is defined for this logical unit.
835 
836  call nemtba(lun,subset,mtyp,msbt,inod)
837 
838  ! Check if there's a long character string in the definition.
839 
840  nte = isc(inod)-inod
841 
842  do i = 1, nte
843  if ( (typ(inod+i)=='CHR') .and. (ibt(inod+i)>64) ) then
844  iret = 1
845  return
846  endif
847  enddo
848 
849  iret = 0
850 
851  return
852 end function lcmgdf
853 
878 recursive subroutine ufbpos(lunit,irec,isub,subset,jdate)
879 
880  use bufrlib
881 
882  use modv_vars, only: im8b
883 
884  use moda_msgcwd
885  use moda_bitbuf
886 
887  implicit none
888 
889  integer, intent(in) :: lunit, irec, isub
890  integer, intent(out) :: jdate
891  integer my_lunit, my_irec, my_isub, lun, il, im, jrec, jsub, iret
892 
893  character*128 bort_str
894  character*8, intent(out) :: subset
895 
896  ! Check for I8 integers
897 
898  if(im8b) then
899  im8b=.false.
900  call x84(lunit,my_lunit,1)
901  call x84(irec,my_irec,1)
902  call x84(isub,my_isub,1)
903  call ufbpos(my_lunit,my_irec,my_isub,subset,jdate)
904  call x48(jdate,jdate,1)
905  im8b=.true.
906  return
907  endif
908 
909  ! Make sure a file is open for input
910 
911  call status(lunit,lun,il,im)
912  if(il==0) call bort('BUFRLIB: UFBPOS - INPUT BUFR FILE IS CLOSED, IT MUST BE OPEN FOR INPUT')
913  if(il>0) call bort('BUFRLIB: UFBPOS - INPUT BUFR FILE IS OPEN FOR OUTPUT, IT MUST BE OPEN FOR INPUT')
914 
915  if(irec<=0) then
916  write(bort_str,'("BUFRLIB: UFBPOS - REQUESTED MESSAGE NUMBER TO READ IN (",I5,") IS NOT VALID")') irec
917  call bort(bort_str)
918  endif
919  if(isub<=0) then
920  write(bort_str,'("BUFRLIB: UFBPOS - REQUESTED SUBSET NUMBER TO READ IN (",I5,") IS NOT VALID")') isub
921  call bort(bort_str)
922  endif
923 
924  ! See where pointers are currently located
925 
926  call ufbcnt(lunit,jrec,jsub)
927 
928  ! Rewind file if requested pointers are behind current pointers
929 
930  if(irec<jrec .or. (irec==jrec.and.isub<jsub)) then
931  call cewind_c(lun)
932  nmsg(lun) = 0
933  nsub(lun) = 0
934  call ufbcnt(lunit,jrec,jsub)
935  endif
936 
937  ! Read subset #isub from message #irec from file
938 
939  do while (irec>jrec)
940  call readmg(lunit,subset,jdate,iret)
941  if(iret<0) then
942  write(bort_str,'("BUFRLIB: UFBPOS - REQUESTED MESSAGE NUMBER '// &
943  'TO READ IN (",I5,") EXCEEDS THE NUMBER OF MESSAGES IN THE FILE (",I5,")")') irec, jrec
944  call bort(bort_str)
945  endif
946  call ufbcnt(lunit,jrec,jsub)
947  enddo
948 
949  do while (isub>jsub)
950  call readsb(lunit,iret)
951  if(iret/=0) then
952  write(bort_str,'("BUFRLIB: UFBPOS - REQ. SUBSET NUMBER TO READ'// &
953  ' IN (",I5,") EXCEEDS THE NUMBER OF SUBSETS (",I5,") IN THE REQ. MESSAGE (",I5,")")') isub, jsub, irec
954  call bort(bort_str)
955  endif
956  call ufbcnt(lunit,jrec,jsub)
957  enddo
958 
959  return
960 end subroutine ufbpos
961 
973 subroutine rdtree(lun,iret)
974 
975  use modv_vars, only: bmiss
976 
977  use moda_usrint
978  use moda_usrbit
979  use moda_ival
980  use moda_bitbuf
981  use moda_tables
982 
983  implicit none
984 
985  integer, intent(in) :: lun
986  integer, intent(out) :: iret
987  integer ier, n, node, kbit, nbt, icbfms, igetrfel
988 
989  character*8 cval
990 
991  real*8 rval, ups
992 
993  equivalence(cval,rval)
994 
995  iret = 0
996 
997  ! Cycle through a subset setting up the template
998 
999  mbit(1) = ibit
1000  nbit(1) = 0
1001  call rcstpl(lun,ier)
1002  if(ier/=0) then
1003  iret = -1
1004  return
1005  endif
1006 
1007  ! Loop through each element of the subset, unpacking each value and then converting it to the proper type
1008 
1009  do n=1,nval(lun)
1010  call upb8(ival(n),nbit(n),mbit(n),mbay(1,lun))
1011  node = inv(n,lun)
1012  if(itp(node)==1) then
1013  ! The unpacked value is a delayed descriptor replication factor.
1014  val(n,lun) = ival(n)
1015  elseif(itp(node)==2) then
1016  ! The unpacked value is a real.
1017  nrfelm(n,lun) = igetrfel(n,lun)
1018  if (ival(n)<2_8**ibt(node)-1) then
1019  val(n,lun) = ups(ival(n),node)
1020  else
1021  val(n,lun) = bmiss
1022  endif
1023  elseif(itp(node)==3) then
1024  ! The value is a character string, so unpack it using an equivalenced real*8 value. Note that a maximum of 8 characters
1025  ! will be unpacked here, so a separate subsequent call to subroutine readlc() will be needed to fully unpack any string
1026  ! longer than 8 characters.
1027  cval = ' '
1028  kbit = mbit(n)
1029  nbt = min(8,nbit(n)/8)
1030  call upc(cval,nbt,mbay(1,lun),kbit,.true.)
1031  if (nbit(n)<=64 .and. icbfms(cval,nbt)/=0) then
1032  val(n,lun) = bmiss
1033  else
1034  val(n,lun) = rval
1035  endif
1036  endif
1037  enddo
1038 
1039  ibit = nbit(nval(lun))+mbit(nval(lun))
1040 
1041  return
1042 end subroutine rdtree
1043 
1052 subroutine wrtree(lun)
1053 
1054  use moda_usrint
1055  use moda_ival
1056  use moda_ufbcpl
1057  use moda_bitbuf
1058  use moda_tables
1059 
1060  implicit none
1061 
1062  integer, intent(in) :: lun
1063  integer*8 ipks
1064  integer n, node, nbit, ncr, numchr, jj, ibfms, igetrfel, imrkopr
1065 
1066  character*120 lstr
1067  character*8 cval
1068 
1069  real*8 rval
1070 
1071  equivalence(cval,rval)
1072 
1073  ! Convert user numbers into scaled integers
1074 
1075  do n=1,nval(lun)
1076  node = inv(n,lun)
1077  nrfelm(n,lun) = igetrfel(n,lun)
1078  if(itp(node)==1) then
1079  ival(n) = nint(val(n,lun))
1080  elseif(typ(node)=='NUM') then
1081  if( (ibfms(val(n,lun))==1) .or. (val(n,lun)/=val(n,lun)) ) then
1082  ! The user number is either "missing" or NaN.
1083  ival(n) = -1
1084  else
1085  ival(n) = ipks(val(n,lun),node)
1086  endif
1087  call strbtm(n,lun,int(ival(n)))
1088  endif
1089  enddo
1090 
1091  ! Pack the user array into the subset buffer
1092 
1093  ibit = 16
1094 
1095  do n=1,nval(lun)
1096  node = inv(n,lun)
1097  if(itp(node)<3) then
1098  ! The value to be packed is numeric.
1099  if ( imrkopr(tag(node)) == 1 ) then
1100  nbit = ibt(inv(nrfelm(n,lun),lun))
1101  else
1102  nbit = ibt(node)
1103  endif
1104  call pkb8(ival(n),nbit,ibay,ibit)
1105  else
1106  ! The value to be packed is a character string.
1107  ncr=ibt(node)/8
1108  if ( ncr>8 .and. luncpy(lun)/=0 ) then
1109  ! The string is longer than 8 characters and there was a preceeding call to ufbcpy() involving this output unit,
1110  ! so read the long string with readlc() and then write it into the output buffer using pkc().
1111  call readlc(luncpy(lun),lstr,tag(node))
1112  call pkc(lstr,ncr,ibay,ibit)
1113  else
1114  rval = val(n,lun)
1115  if(ibfms(rval)/=0) then
1116  ! The value is "missing", so set all bits to 1 before packing the field as a character string.
1117  numchr = min(ncr,len(lstr))
1118  do jj = 1, numchr
1119  call ipkm(lstr(jj:jj),1,255)
1120  enddo
1121  call pkc(lstr,numchr,ibay,ibit)
1122  else
1123  ! The value is not "missing", so pack the equivalenced character string. Note that a maximum of 8 characters
1124  ! will be packed here, so a separate subsequent call to subroutine writlc() will be needed to fully encode any
1125  ! string longer than 8 characters.
1126  call pkc(cval,ncr,ibay,ibit)
1127  endif
1128  endif
1129  endif
1130  enddo
1131 
1132  ! Reset ufbcpy() file pointer
1133 
1134  luncpy(lun)=0
1135 
1136  return
1137 end subroutine wrtree
1138 
1151 subroutine rcstpl(lun,iret)
1152 
1153  use modv_vars, only: maxjl, maxss, maxrcr, iprt
1154 
1155  use moda_usrint
1156  use moda_usrbit
1157  use moda_msgcwd
1158  use moda_bitbuf
1159  use moda_tables
1160  use moda_usrtmp
1161 
1162  implicit none
1163 
1164  character*128 bort_str
1165 
1166  integer, intent(in) :: lun
1167  integer, intent(out) :: iret
1168  integer nbmp(2,maxrcr), newn(2,maxrcr), knx(maxrcr), nodi, node, mbmp, nr, i, j, n, nn, n1, n2, new, ivob, igetrfel
1169 
1170  iret = 0
1171 
1172  ! Set the initial values for the template
1173 
1174  inv(1,lun) = inode(lun)
1175  val(1,lun) = 0
1176  nbmp(1,1) = 1
1177  nbmp(2,1) = 1
1178  nodi = inode(lun)
1179  node = inode(lun)
1180  mbmp = 1
1181  nval(lun) = 1
1182  nr = 0
1183  knx(1:maxrcr) = 0
1184 
1185  outer: do while (.true.)
1186 
1187  ! Set up the parameters for a level of recursion
1188 
1189  nr = nr+1
1190  if(nr>maxrcr) then
1191  write(bort_str,'("BUFRLIB: RCSTPL - THE NUMBER OF RECURSION LEVELS EXCEEDS THE LIMIT (",I3,")")') maxrcr
1192  call bort(bort_str)
1193  endif
1194  nbmp(1,nr) = 1
1195  nbmp(2,nr) = mbmp
1196 
1197  n1 = iseq(node,1)
1198  n2 = iseq(node,2)
1199  if(n1==0) then
1200  write(bort_str,'("BUFRLIB: RCSTPL - UNSET EXPANSION SEGMENT ",A)') tag(nodi)
1201  call bort(bort_str)
1202  endif
1203  if(n2-n1+1>maxjl) then
1204  if(iprt>=0) then
1205  call errwrt('++++++++++++++BUFR ARCHIVE LIBRARY+++++++++++++++++')
1206  call errwrt('BUFRLIB: RCSTPL - MAXJL OVERFLOW; SUBSET SKIPPED')
1207  call errwrt('++++++++++++++BUFR ARCHIVE LIBRARY+++++++++++++++++')
1208  endif
1209  iret = -1
1210  return
1211  endif
1212  newn(1,nr) = 1
1213  newn(2,nr) = n2-n1+1
1214 
1215  do n=1,newn(2,nr)
1216  nn = jseq(n+n1-1)
1217  iutmp(n,nr) = nn
1218  vutmp(n,nr) = vali(nn)
1219  enddo
1220 
1221  do while (.true.)
1222 
1223  ! Store nodes at some recursion level
1224 
1225  do i=nbmp(1,nr),nbmp(2,nr)
1226  if(knx(nr)==0) knx(nr) = nval(lun)
1227  if(i>nbmp(1,nr)) newn(1,nr) = 1
1228  do j=newn(1,nr),newn(2,nr)
1229  if(nval(lun)+1>maxss) then
1230  if(iprt>=0) then
1231  call errwrt('++++++++++++++BUFR ARCHIVE LIBRARY+++++++++++++++++')
1232  call errwrt('BUFRLIB: RCSTPL - MAXSS OVERFLOW; SUBSET SKIPPED')
1233  call errwrt('++++++++++++++BUFR ARCHIVE LIBRARY+++++++++++++++++')
1234  endif
1235  iret = -1
1236  return
1237  endif
1238  nval(lun) = nval(lun)+1
1239  node = iutmp(j,nr)
1240  ! inv is positional index in internal jump/link table for packed subset element nval(lun) in mbay
1241  inv(nval(lun),lun) = node
1242  ! mbit is the bit in mbay pointing to where the packed subset element nval(lun) begins
1243  mbit(nval(lun)) = mbit(nval(lun)-1)+nbit(nval(lun)-1)
1244  ! nbit is the number of bits in mbay occupied by packed subset element nval(lun)
1245  nrfelm(nval(lun),lun) = igetrfel(nval(lun),lun)
1246  nbit(nval(lun)) = ibt(node)
1247  if(nbit(nval(lun))==1) then
1248  ! Check whether this is a bitmap entry
1249  call upbb(ivob,nbit(nval(lun)),mbit(nval(lun)),mbay(1,lun))
1250  call strbtm(nval(lun),lun,ivob)
1251  endif
1252  ! Actual unpacked subset values are initialized here
1253  val(nval(lun),lun) = vutmp(j,nr)
1254  if(itp(node)==1) then
1255  call upbb(mbmp,nbit(nval(lun)),mbit(nval(lun)),mbay(1,lun))
1256  newn(1,nr) = j+1
1257  nbmp(1,nr) = i
1258  cycle outer
1259  endif
1260  enddo
1261  new = nval(lun)-knx(nr)
1262  val(knx(nr)+1,lun) = val(knx(nr)+1,lun) + new
1263  knx(nr) = 0
1264  enddo
1265 
1266  ! Check if we need to continue one recursion level back
1267 
1268  if(nr-1 == 0) exit outer
1269  nr = nr-1
1270  enddo
1271 
1272  enddo outer
1273 
1274  return
1275 end subroutine rcstpl
1276 
1287 subroutine usrtpl(lun,invn,nbmp)
1288 
1289  use modv_vars, only: maxjl, maxss, iprt
1290 
1291  use moda_usrint
1292  use moda_msgcwd
1293  use moda_tables
1294  use moda_ivttmp
1295  use moda_stcode
1296 
1297  implicit none
1298 
1299  integer, intent(in) :: lun, invn, nbmp
1300  integer i, j, ival, jval, n, n1, n2, nodi, node, newn, invr, knvn
1301 
1302  character*128 bort_str, errstr
1303 
1304  logical drp, drs, drb, drx
1305 
1306  if(iprt>=2) then
1307  call errwrt('++++++++++++++BUFR ARCHIVE LIBRARY+++++++++++++++++')
1308  write ( unit=errstr, fmt='(A,I3,A,I7,A,I5,A,A10)' ) &
1309  'BUFRLIB: USRTPL - LUN:INVN:NBMP:TAG(INODE(LUN)) = ', lun, ':', invn, ':', nbmp, ':', tag(inode(lun))
1310  call errwrt(errstr)
1311  call errwrt('++++++++++++++BUFR ARCHIVE LIBRARY+++++++++++++++++')
1312  call errwrt(' ')
1313  endif
1314 
1315  if(nbmp<=0) then
1316  if(iprt>=1) then
1317  call errwrt('+++++++++++++++++++++WARNING+++++++++++++++++++++++')
1318  call errwrt(.LE.'BUFRLIB: USRTPL - NBMP 0 - IMMEDIATE RETURN')
1319  call errwrt('+++++++++++++++++++++WARNING+++++++++++++++++++++++')
1320  call errwrt(' ')
1321  endif
1322  return
1323  endif
1324 
1325  drp = .false.
1326  drs = .false.
1327  drx = .false.
1328 
1329  ! Set up a node expansion
1330 
1331  if(invn==1) then
1332  ! The node is a Table A mnemonic
1333  nodi = inode(lun)
1334  inv(1,lun) = nodi
1335  nval(lun) = 1
1336  if(nbmp/=1) then
1337  write(bort_str,'("BUFRLIB: USRTPL - THIRD ARGUMENT (INPUT) = ",'// &
1338  'I4,", MUST BE 1 WHEN SECOND ARGUMENT (INPUT) IS 1 (SUBSET NODE) (",A,")")') nbmp, tag(nodi)
1339  call bort(bort_str)
1340  endif
1341  elseif(invn>0 .and. invn<=nval(lun)) then
1342  ! The node is (hopefully) a delayed replication factor
1343  nodi = inv(invn,lun)
1344  drp = typ(nodi) == 'DRP'
1345  drs = typ(nodi) == 'DRS'
1346  drb = typ(nodi) == 'DRB'
1347  drx = drp .or. drs .or. drb
1348  ival = nint(val(invn,lun))
1349  jval = 2**ibt(nodi)-1
1350  val(invn,lun) = ival+nbmp
1351  if(drb.and.nbmp/=1) then
1352  write(bort_str,'("BUFRLIB: USRTPL - THIRD ARGUMENT (INPUT) = ",'// &
1353  'I4,", MUST BE 1 WHEN NODE IS DRB (1-BIT DELAYED REPL. FACTOR) (",A,")")') nbmp, tag(nodi)
1354  call bort(bort_str)
1355  endif
1356  if(.not.drx) then
1357  write(bort_str,'("BUFRLIB: USRTPL - NODE IS OF TYPE ",A," - IT '// &
1358  'MUST BE EITHER A SUBSET OR DELAYED REPL. FACTOR (",A,")")') typ(nodi), tag(nodi)
1359  call bort(bort_str)
1360  endif
1361  if(ival<0) then
1362  write(bort_str,'("BUFRLIB: USRTPL - REPLICATION FACTOR IS NEGATIVE (=",I5,") (",A,")")') ival, tag(nodi)
1363  call bort(bort_str)
1364  endif
1365  if(ival+nbmp>jval) then
1366  write(bort_str,'("BUFRLIB: USRTPL - REPLICATION FACTOR OVERFLOW (EXCEEDS MAXIMUM OF",I6," (",A,")")') jval, tag(nodi)
1367  call errwrt(bort_str)
1368  iscodes(lun) = 1
1369  return
1370  endif
1371  else
1372  write(bort_str,'("BUFRLIB: USRTPL - INVENTORY INDEX {FIRST '// &
1373  'ARGUMENT (INPUT)} OUT OF BOUNDS (=",I5,", RANGE IS 1 TO",I6,") ")') invn, nval(lun)
1374  call bort(bort_str)
1375  endif
1376 
1377  ! Recall a pre-fab node expansion segment
1378 
1379  newn = 0
1380  n1 = iseq(nodi,1)
1381  n2 = iseq(nodi,2)
1382 
1383  if(n1==0) then
1384  write(bort_str,'("BUFRLIB: USRTPL - UNSET EXPANSION SEGMENT (",A,")")') tag(nodi)
1385  call bort(bort_str)
1386  endif
1387  if(n2-n1+1>maxjl) then
1388  write(bort_str,'("BUFRLIB: USRTPL - TEMPLATE ARRAY OVERFLOW, EXCEEDS THE LIMIT (",I6,") (",A,")")') maxjl, tag(nodi)
1389  call bort(bort_str)
1390  endif
1391 
1392  do n=n1,n2
1393  newn = newn+1
1394  itmp(newn) = jseq(n)
1395  vtmp(newn) = vali(jseq(n))
1396  enddo
1397 
1398  ! Move old nodes and store new ones
1399 
1400  if(nval(lun)+newn*nbmp>maxss) then
1401  write(bort_str,'("BUFRLIB: USRTPL - INVENTORY OVERFLOW (",I6,"), EXCEEDS THE LIMIT (",I6,") (",A,")")') &
1402  nval(lun)+newn*nbmp, maxss, tag(nodi)
1403  call bort(bort_str)
1404  endif
1405 
1406  do j=nval(lun),invn+1,-1
1407  inv(j+newn*nbmp,lun) = inv(j,lun)
1408  val(j+newn*nbmp,lun) = val(j,lun)
1409  enddo
1410 
1411  if(drp.or.drs) vtmp(1) = newn
1412  knvn = invn
1413 
1414  do i=1,nbmp
1415  do j=1,newn
1416  knvn = knvn+1
1417  inv(knvn,lun) = itmp(j)
1418  val(knvn,lun) = vtmp(j)
1419  enddo
1420  enddo
1421 
1422  ! Reset pointers and counters
1423 
1424  nval(lun) = nval(lun) + newn*nbmp
1425 
1426  if(iprt>=2) then
1427  call errwrt('++++++++++++++BUFR ARCHIVE LIBRARY+++++++++++++++++')
1428  write ( unit=errstr, fmt='(A,A,A10,2(A,I5),A,I7)' ) 'BUFRLIB: USRTPL - TAG(INV(INVN,LUN)):NEWN:NBMP:', &
1429  'NVAL(LUN) = ', tag(inv(invn,lun)), ':', newn, ':', nbmp, ':', nval(lun)
1430  call errwrt(errstr)
1431  do i=1,newn
1432  write ( unit=errstr, fmt='(2(A,I5),A,A10)' ) 'For I = ', i, ', ITMP(I) = ', itmp(i), ', TAG(ITMP(I)) = ', tag(itmp(i))
1433  call errwrt(errstr)
1434  enddo
1435  call errwrt('++++++++++++++BUFR ARCHIVE LIBRARY+++++++++++++++++')
1436  call errwrt(' ')
1437  endif
1438 
1439  if(drx) then
1440  node = nodi
1441  invr = invn
1442  outer: do while (.true.)
1443  node = jmpb(node)
1444  if(node<=0) exit
1445  if(itp(node)==0) then
1446  do invr=invr-1,1,-1
1447  if(inv(invr,lun)==node) then
1448  val(invr,lun) = val(invr,lun)+newn*nbmp
1449  cycle outer
1450  endif
1451  enddo
1452  write(bort_str,'("BUFRLIB: USRTPL - BAD BACKUP STRATEGY (",A,")")') tag(nodi)
1453  call bort(bort_str)
1454  else
1455  cycle
1456  endif
1457  enddo outer
1458  endif
1459 
1460  return
1461 end subroutine usrtpl
1462 
1476 recursive subroutine invmrg(lubfi,lubfj)
1477 
1478  use modv_vars, only: im8b
1479 
1480  use moda_usrint
1481  use moda_tables
1482  use moda_mrgcom
1483 
1484  implicit none
1485 
1486  integer, intent(in) :: lubfi, lubfj
1487  integer my_lubfi, my_lubfj, luni, il, im, lunj, jl, jm, is, js, node, nodj, ityp, iwrds, jwrds, &
1488  n, ioff, nwords, ibfms
1489 
1490  character*128 bort_str
1491 
1492  logical herei, herej, missi, missj, samei
1493 
1494  ! Check for I8 integers
1495 
1496  if(im8b) then
1497  im8b=.false.
1498  call x84(lubfi,my_lubfi,1)
1499  call x84(lubfj,my_lubfj,1)
1500  call invmrg(my_lubfi,my_lubfj)
1501  im8b=.true.
1502  return
1503  endif
1504 
1505  is = 1
1506  js = 1
1507 
1508  ! Get the unit pointers
1509 
1510  call status(lubfi,luni,il,im)
1511  call status(lubfj,lunj,jl,jm)
1512 
1513  ! Step through the buffers comparing the inventory and merging data
1514 
1515  do while(is<=nval(luni))
1516  ! Confirm we're at the same node in each buffer
1517  node = inv(is,luni)
1518  nodj = inv(js,lunj)
1519  if(node/=nodj) then
1520  write(bort_str,'("BUFRLIB: INVMRG - NODE FROM INPUT BUFR FILE '// &
1521  '(",I7,") DOES NOT EQUAL NODE FROM OUTPUT BUFR FILE (",I7,"), TABULAR MISMATCH")') node, nodj
1522  call bort(bort_str)
1523  endif
1524 
1525  ityp = itp(node)
1526  if(ityp==1) then
1527  ! Do an entire sequence replacement
1528  if(typ(node)=='DRB') then
1529  ioff = 0
1530  else
1531  ioff = 1
1532  endif
1533  iwrds = nwords(is,luni)+ioff
1534  jwrds = nwords(js,lunj)+ioff
1535  if(iwrds>ioff .and. jwrds==ioff) then
1536  do n=nval(lunj),js+1,-1
1537  inv(n+iwrds-jwrds,lunj) = inv(n,lunj)
1538  val(n+iwrds-jwrds,lunj) = val(n,lunj)
1539  enddo
1540  do n=0,iwrds
1541  inv(js+n,lunj) = inv(is+n,luni)
1542  val(js+n,lunj) = val(is+n,luni)
1543  enddo
1544  nval(lunj) = nval(lunj)+iwrds-jwrds
1545  jwrds = iwrds
1546  nrpl = nrpl+1
1547  endif
1548  is = is+iwrds
1549  js = js+jwrds
1550  elseif((ityp==2).or.(ityp==3)) then
1551  ! Fill missing values
1552  herei = ibfms(val(is,luni))==0
1553  herej = ibfms(val(js,lunj))==0
1554  missi = .not.(herei)
1555  missj = .not.(herej)
1556  samei = val(is,luni)==val(js,lunj)
1557  if(herei.and.missj) then
1558  val(js,lunj) = val(is,luni)
1559  nmrg = nmrg+1
1560  elseif(herei.and.herej.and..not.samei) then
1561  namb = namb+1
1562  endif
1563  endif
1564 
1565  ! Bump the counters and go check the next pair
1566  is = is + 1
1567  js = js + 1
1568  enddo
1569 
1570  ntot = ntot+1
1571 
1572  return
1573 end subroutine invmrg
1574 
1583 integer function nwords(n,lun) result(iret)
1584 
1585  use moda_usrint
1586 
1587  implicit none
1588 
1589  integer, intent(in) :: n, lun
1590  integer k
1591 
1592  iret = 0
1593 
1594  do k=1,nint(val(n,lun))
1595  iret = iret + nint(val(iret+n+1,lun))
1596  enddo
1597 
1598  return
1599 end function nwords
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
subroutine bort(str)
Log an error message, then abort the application program.
Definition: borts.F90:15
subroutine upb(nval, nbits, ibay, ibit)
Decode an integer value from within a specified number of bits of an integer array,...
Definition: cidecode.F90:202
subroutine upbb(nval, nbits, ibit, ibay)
Decode an integer value from within a specified number of bits of an integer array,...
Definition: cidecode.F90:154
subroutine upb8(nval, nbits, ibit, ibay)
Decode an 8-byte integer value from within a specified number of bits of an integer array,...
Definition: cidecode.F90:80
real *8 function ups(ival, node)
Unpack a real*8 value from an integer by applying the proper scale and reference values.
Definition: cidecode.F90:319
subroutine upc(chr, nchr, ibay, ibit, cnvnull)
Decode a character string from within a specified number of bytes of an integer array,...
Definition: cidecode.F90:26
subroutine pkc(chr, nchr, ibay, ibit)
Encode a character string within a specified number of bytes of an integer array, starting at the bit...
Definition: ciencode.F90:25
recursive subroutine ipkm(cbay, nbyt, n)
Encode an integer value within a specified number of bytes of a character string, up to a maximum of ...
Definition: ciencode.F90:194
subroutine pkb(nval, nbits, ibay, ibit)
Encode an integer value within a specified number of bits of an integer array, starting at the bit im...
Definition: ciencode.F90:140
subroutine pkb8(nval, nbits, ibay, ibit)
Encode an 8-byte integer value within a specified number of bits of an integer array,...
Definition: ciencode.F90:97
subroutine rdcmps(lun)
Read the next compressed BUFR data subset into internal arrays.
Definition: compress.F90:110
subroutine wrcmps(lunix)
Write a compressed BUFR data subset.
Definition: compress.F90:384
subroutine mvb(ib1, nb1, ib2, nb2, nbm)
Copy a specified number of bytes from one packed binary array to another.
Definition: copydata.F90:729
subroutine nemtba(lun, nemo, mtyp, msbt, inod)
Get information about a Table A descriptor from the internal DX BUFR tables.
Definition: dxtable.F90:1238
subroutine errwrt(str)
Specify a custom location for the logging of error and diagnostic messages generated by the NCEPLIBS-...
Definition: errwrt.F90:32
integer function ibfms(r8val)
Check whether a real*8 data value returned from a previous call to any of the NCEPLIBS-bufr values-re...
Definition: missing.F90:25
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 BUFR messages internally for multiple file IDs.
integer, dimension(:), allocatable ibay
Current data subset.
integer ibit
Bit pointer within ibay.
integer, dimension(:,:), allocatable mbay
Current BUFR message for each file ID.
integer, dimension(:), allocatable mbyt
Length (in bytes) of current BUFR message for each file ID.
integer maxbyt
Maximum length of an output BUFR message.
Declare arrays and variables used to store bitmaps internally within a data subset definition.
integer lstnod
Most recent jump/link table entry that was processed by function igetrfel() and whose corresponding v...
integer nbtm
Number of stored bitmaps for the current data subset (up to a maximum of mxbtm).
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.
Declare arrays used to store, for each output file ID, a copy of the BUFR message that was most recen...
integer, dimension(:), allocatable msglen
Length (in integers) of BUFR message most recently written to each output file ID.
integer, dimension(:,:), allocatable msgtxt
BUFR message most recently written to each output file ID.
Declare arrays and variables needed to store long character strings (greater than 8 bytes) via subrou...
integer nh4wlc
Number of long character strings being stored.
character *14, dimension(:), allocatable sth4wlc
Table B mnemonics associated with long character strings.
integer, dimension(:), allocatable luh4wlc
File ID for associated output file.
character *120, dimension(:), allocatable chh4wlc
Long character strings.
Declare an array used to pack or unpack all of the values of a BUFR data subset.
integer *8, dimension(:), allocatable ival
BUFR data subset values.
Declare arrays which provide working space in several subprograms (usrtpl() and ufbcup()) which manip...
real *8, dimension(:), allocatable vtmp
val array elements for new sections of a growing subset buffer.
integer, dimension(:), allocatable itmp
inv array elements for new sections of a growing subset buffer.
Declare variables for use when merging parts of different data subsets.
integer nmrg
Number of merges.
integer ntot
Total number of calls to subroutine invmrg().
integer namb
Number of potential merges that weren't made because of ambiguities.
integer nrpl
Number of expansions of Table D mnemonics using short (1-bit) delayed replication.
Declare a variable used to indicate whether output BUFR messages should be compressed.
character ccmf
Flag indicating whether BUFR output messages are to be compressed; this variable is initialized to a ...
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.
integer, dimension(:), allocatable idate
Section 1 date-time of message.
integer, dimension(:), allocatable nmsg
Current message pointer within logical unit.
integer, dimension(:), allocatable msub
Total number of data subsets in message.
integer, dimension(:), allocatable nsub
Current subset pointer within message.
Declare an array used to store a status code for each file ID if an error or other abnormal result oc...
integer, dimension(:), allocatable iscodes
Abnormal status codes.
Declare arrays and variables used to store the internal jump/link table.
integer, dimension(:), allocatable jseq
Temporary storage used in expanding sequences.
integer, dimension(:,:), allocatable iseq
Temporary storage used in expanding sequences.
integer, dimension(:), allocatable isc
Scale factors corresponding to tag and typ:
integer, dimension(:), allocatable ibt
Bit widths corresponding to tag and typ:
real *8, dimension(:), allocatable vali
Initialized data values corresponding to typ:
character *3, dimension(:), allocatable typ
Type indicators corresponding to tag:
integer, dimension(:), allocatable jmpb
Jump backward indices 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 an array used to store, for each file ID, the logical unit number corresponding to a separate...
integer, dimension(:), allocatable luncpy
Logical unit numbers used to copy long character strings between BUFR data subsets.
Declare an array used to store, for each file ID from which a BUFR message is currently being read as...
integer, dimension(:), allocatable msgunp
Flag indicating how to unpack data subsets from BUFR message:
Declare arrays for internal storage of pointers to BUFR data subset values.
integer, dimension(:), allocatable nbit
Length (in bits) of each packed data value in data subset.
integer, dimension(:), allocatable mbit
Pointer in data subset to first bit of each packed data value.
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...
Declare arrays used in subroutine rcstpl() to store subset segments that are being copied from a subs...
integer, dimension(:,:), allocatable iutmp
inv array elements for new sections of a growing subset buffer.
real *8, dimension(:,:), allocatable vutmp
val array elements for new sections of a growing subset buffer.
recursive subroutine openbf(lunit, io, lundx)
Connect a new file to the NCEPLIBS-bufr software for input or output operations, or initialize the li...
recursive subroutine ufbcnt(lunit, kmsg, ksub)
Get the current location of the file pointer within a BUFR file, in terms of a message number countin...
recursive subroutine status(lunit, lun, il, im)
Check whether a specified Fortran logical unit number is currently connected to the NCEPLIBS-bufr sof...
recursive subroutine closmg(lunin)
Close the BUFR message that is currently open for writing within internal arrays associated with logi...
recursive subroutine readmg(lunxx, subset, jdate, iret)
Read the next BUFR message from logical unit abs(lunxx) into internal arrays.
Definition: readwritemg.F90:44
subroutine msgini(lun)
Initialize, within the internal arrays, a new uncompressed BUFR message for output.
subroutine msgwrt(lunit, mesg, mgbyt)
Perform final checks and updates on a BUFR message before writing it to a specified Fortran logical u...
recursive subroutine ufbpos(lunit, irec, isub, subset, jdate)
Jump forwards or backwards to a specified data subset within a BUFR file.
subroutine pad(ibay, ibit, ibyt, ipadb)
Pad a BUFR data subset with zeroed-out bits up to the next byte boundary.
subroutine rdtree(lun, iret)
Read the next uncompressed BUFR data subset into internal arrays.
subroutine wrtree(lun)
Pack a BUFR data subset.
subroutine msgupd(lunit, lun)
Write an uncompressed BUFR data subset.
subroutine usrtpl(lun, invn, nbmp)
Expand a subset template within internal arrays.
recursive subroutine writsa(lunxx, lmsgt, msgt, msgl)
Write a complete data subset into a BUFR message, and return each completed message within a memory a...
recursive integer function ireadns(lunit, subset, idate)
Call subroutine readns() and pass back its return code as the function value.
integer function nwords(n, lun)
Compute the length of a specified delayed replication sequence within a data subset.
recursive subroutine writsb(lunit)
Write a complete data subset into a BUFR message, for eventual output to logical unit lunit.
recursive subroutine invmrg(lubfi, lubfj)
Merge parts of data subsets which have duplicate space and time coordinates but different or unique o...
recursive integer function ireadsb(lunit)
Call subroutine readsb() and pass back its return code as the function value.
recursive subroutine readsb(lunit, iret)
Read the next data subset from a BUFR message.
Definition: readwritesb.F90:32
subroutine rcstpl(lun, iret)
Initialize a subset template within internal arrays.
recursive integer function lcmgdf(lunit, subset)
Check whether the subset definition for a given message type contains any long character strings (gre...
recursive subroutine rdmgsb(lunit, imsg, isub)
Read a specified data subset from a BUFR file.
recursive subroutine readns(lunit, subset, jdate, iret)
Read the next data subset from a BUFR file.
recursive subroutine readlc(lunit, chr, str)
Read a long character string (greater than 8 bytes) from a data subset.
recursive subroutine writlc(lunit, chr, str)
Write a long character string (greater than 8 bytes) to a data subset.
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