NCEPLIBS-g2  3.4.8
getgb2rp.F90
Go to the documentation of this file.
1 
5 
38 subroutine getgb2rp(lugb, cindex, extract, gribm, leng, iret)
39  implicit none
40 
41  integer, intent(in) :: lugb
42  character(len = 1), intent(in) :: cindex(*)
43  logical, intent(in) :: extract
44  integer, intent(out) :: leng, iret
45  character(len = 1), pointer, dimension(:) :: gribm
46 
47  integer, parameter :: zero = 0
48  character(len = 1), allocatable, dimension(:) :: csec2, csec6, csec7
49  character(len = 4) :: ctemp
50  integer :: lencur, lread, len0, ibmap, ipos, iskip
51  integer :: len7, len8, len3, len4, len5, len6, len1, len2
52  integer :: iskp2, iskp6, iskp7
53 
54  iret = 0
55 
56  ! Extract grib message from file.
57  IF (extract) THEN
58  len0 = 16
59  len8 = 4
60  CALL g2_gbytec(cindex, iskip, 4*8, 4*8) ! BYTES TO SKIP IN FILE
61  CALL g2_gbytec(cindex, iskp2, 8*8, 4*8) ! BYTES TO SKIP FOR section 2
62  if (iskp2 .gt. 0) then
63  CALL baread(lugb, iskip + iskp2, 4, lread, ctemp)
64  CALL g2_gbytec(ctemp, len2, 0, 4*8) ! LENGTH OF SECTION 2
65  ALLOCATE(csec2(len2))
66  CALL baread(lugb, iskip + iskp2, len2, lread, csec2)
67  else
68  len2 = 0
69  endif
70  CALL g2_gbytec(cindex, len1, 44*8, 4*8) ! LENGTH OF SECTION 1
71  ipos = 44 + len1
72  CALL g2_gbytec(cindex, len3, ipos*8, 4*8) ! LENGTH OF SECTION 3
73  ipos = ipos + len3
74  CALL g2_gbytec(cindex, len4, ipos*8, 4*8) ! LENGTH OF SECTION 4
75  ipos = ipos + len4
76  CALL g2_gbytec(cindex, len5, ipos*8, 4*8) ! LENGTH OF SECTION 5
77  ipos = ipos + len5
78  CALL g2_gbytec(cindex, len6, ipos*8, 4*8) ! LENGTH OF SECTION 6
79  ipos = ipos + 5
80  CALL g2_gbytec(cindex, ibmap, ipos*8, 1*8) ! Bitmap indicator
81  IF (ibmap .eq. 254) THEN
82  CALL g2_gbytec(cindex, iskp6, 24*8, 4*8) ! BYTES TO SKIP FOR section 6
83  CALL baread(lugb, iskip + iskp6, 4, lread, ctemp)
84  CALL g2_gbytec(ctemp, len6, 0, 4*8) ! LENGTH OF SECTION 6
85  ENDIF
86 
87  ! READ IN SECTION 7 from file
88  CALL g2_gbytec(cindex, iskp7, 28*8, 4*8) ! BYTES TO SKIP FOR section 7
89  CALL baread(lugb, iskip + iskp7, 4, lread, ctemp)
90  CALL g2_gbytec(ctemp, len7, 0, 4*8) ! LENGTH OF SECTION 7
91  ALLOCATE(csec7(len7))
92  CALL baread(lugb, iskip + iskp7, len7, lread, csec7)
93 
94  leng = len0 + len1 + len2 + len3 + len4 + len5 + len6 + len7 + len8
95  IF (.NOT. ASSOCIATED(gribm)) ALLOCATE(gribm(leng))
96 
97  ! Create Section 0
98  gribm(1) = 'G'
99  gribm(2) = 'R'
100  gribm(3) = 'I'
101  gribm(4) = 'B'
102  gribm(5) = char(0)
103  gribm(6) = char(0)
104  gribm(7) = cindex(42)
105  gribm(8) = cindex(41)
106  gribm(9) = char(0)
107  gribm(10) = char(0)
108  gribm(11) = char(0)
109  gribm(12) = char(0)
110  CALL g2_sbytec(gribm, leng, 12*8, 4*8)
111 
112  ! Copy Section 1
113  gribm(17:16 + len1) = cindex(45:44 + len1)
114  lencur = 16 + len1
115  ipos = 44 + len1
116 
117  ! Copy Section 2, if necessary
118  if (iskp2 .gt. 0) then
119  gribm(lencur + 1:lencur + len2) = csec2(1:len2)
120  lencur = lencur + len2
121  endif
122 
123  ! Copy Sections 3 through 5
124  gribm(lencur + 1:lencur + len3 + len4 + len5) = cindex(ipos + 1:ipos + len3 + len4 + len5)
125  lencur = lencur + len3 + len4 + len5
126  ipos = ipos + len3 + len4 + len5
127 
128  ! Copy Section 6
129  if (len6 .eq. 6 .AND. ibmap .ne. 254) then
130  gribm(lencur + 1:lencur + len6) = cindex(ipos + 1:ipos + len6)
131  lencur = lencur + len6
132  else
133  CALL g2_gbytec(cindex, iskp6, 24*8, 4*8) ! BYTES TO SKIP FOR section 6
134  CALL baread(lugb, iskip + iskp6, 4, lread, ctemp)
135  CALL g2_gbytec(ctemp, len6, 0, 4*8) ! LENGTH OF SECTION 6
136  ALLOCATE(csec6(len6))
137  CALL baread(lugb, iskip + iskp6, len6, lread, csec6)
138  gribm(lencur + 1:lencur + len6) = csec6(1:len6)
139  lencur = lencur + len6
140  IF (allocated(csec6)) DEALLOCATE(csec6)
141  endif
142 
143  ! Copy Section 7
144  gribm(lencur + 1:lencur + len7) = csec7(1:len7)
145  lencur = lencur + len7
146 
147  ! Section 8
148  gribm(lencur + 1) = '7'
149  gribm(lencur + 2) = '7'
150  gribm(lencur + 3) = '7'
151  gribm(lencur + 4) = '7'
152 
153  ! clean up
154  IF (allocated(csec2)) DEALLOCATE(csec2)
155  IF (allocated(csec7)) deallocate(csec7)
156  ELSE ! DO NOT extract field from message : Get entire message
157  CALL g2_gbytec(cindex, iskip, 4*8, 4*8) ! BYTES TO SKIP IN FILE
158  CALL g2_gbytec(cindex, leng, 36*8, 4*8) ! LENGTH OF GRIB MESSAGE
159  IF (.NOT. ASSOCIATED(gribm)) ALLOCATE(gribm(leng))
160  CALL baread(lugb, iskip, leng, lread, gribm)
161  IF (leng .NE. lread ) THEN
162  DEALLOCATE(gribm)
163  NULLIFY(gribm)
164  iret = 97
165  RETURN
166  ENDIF
167  ENDIF
168 END SUBROUTINE getgb2rp
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
subroutine getgb2rp(lugb, cindex, extract, gribm, leng, iret)
Extract a grib message from a file given the index of the requested field.
Definition: getgb2rp.F90:39