NCEPLIBS-g2  3.4.7
gf_unpack5.F90
Go to the documentation of this file.
1 
7 
33 subroutine gf_unpack5(cgrib, lcgrib, iofst, ndpts, idrsnum, idrstmpl, &
34  mapdrslen, ierr)
35  use drstemplates
36  use re_alloc ! needed for subroutine realloc
37  implicit none
38 
39  character(len = 1), intent(in) :: cgrib(lcgrib)
40  integer, intent(in) :: lcgrib
41  integer, intent(inout) :: iofst
42  integer, intent(out) :: ndpts, idrsnum
43  integer, pointer, dimension(:) :: idrstmpl
44  integer, intent(out) :: ierr
45 
46  integer, allocatable :: mapdrs(:)
47  integer :: mapdrslen
48  logical needext
49  integer :: newmapdrslen, nbits, istat, isign, lensec, iret, i
50 
51  ierr = 0
52  nullify(idrstmpl)
53 
54  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
55  iofst = iofst + 32
56  iofst = iofst + 8 ! skip section number
57  allocate(mapdrs(lensec))
58 
59  ! Get num of data points.
60  call g2_gbytec(cgrib, ndpts, iofst, 32)
61  iofst = iofst + 32
62  ! Get Data Rep Template Num.
63  call g2_gbytec(cgrib, idrsnum, iofst, 16)
64  iofst = iofst + 16
65  ! Gen Data Representation Template.
66  call getdrstemplate(idrsnum, mapdrslen, mapdrs, needext, iret)
67  if (iret .ne. 0) then
68  ierr = 7
69  if (allocated(mapdrs)) deallocate(mapdrs)
70  return
71  endif
72 
73  ! Unpack each value into array ipdstmpl from the the appropriate
74  ! number of octets, which are specified in corresponding entries in
75  ! array mappds.
76  istat = 0
77  if (mapdrslen .gt. 0) allocate(idrstmpl(mapdrslen), stat = istat)
78  if (istat .ne. 0) then
79  ierr = 6
80  nullify(idrstmpl)
81  if (allocated(mapdrs)) deallocate(mapdrs)
82  return
83  endif
84  do i = 1, mapdrslen
85  nbits = iabs(mapdrs(i)) * 8
86  if (mapdrs(i) .ge. 0) then
87  call g2_gbytec(cgrib, idrstmpl(i), iofst, nbits)
88  else
89  call g2_gbytec(cgrib, isign, iofst, 1)
90  call g2_gbytec(cgrib, idrstmpl(i), iofst + 1, nbits-1)
91  if (isign .eq. 1) idrstmpl(i) = -idrstmpl(i)
92  endif
93  iofst = iofst + nbits
94  enddo
95 
96  ! Check to see if the Data Representation Template needs to be
97  ! extended. The number of values in a specific template may vary
98  ! depending on data specified in the "static" part of the template.
99  if (needext) then
100  call extdrstemplate(idrsnum, idrstmpl, newmapdrslen, mapdrs)
101  call realloc(idrstmpl, mapdrslen, newmapdrslen, istat)
102 
103  ! Unpack the rest of the Data Representation Template.
104  do i = mapdrslen + 1, newmapdrslen
105  nbits = iabs(mapdrs(i)) * 8
106  if (mapdrs(i) .ge. 0) then
107  call g2_gbytec(cgrib, idrstmpl(i), iofst, nbits)
108  else
109  call g2_gbytec(cgrib, isign, iofst, 1)
110  call g2_gbytec(cgrib, idrstmpl(i), iofst + 1, nbits - 1)
111  if (isign.eq.1) idrstmpl(i) = -idrstmpl(i)
112  endif
113  iofst = iofst + nbits
114  enddo
115  mapdrslen = newmapdrslen
116  endif
117  if (allocated(mapdrs)) deallocate(mapdrs)
118 
119 end subroutine gf_unpack5
120 
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 gf_unpack5(cgrib, lcgrib, iofst, ndpts, idrsnum, idrstmpl, mapdrslen, ierr)
Unpack Section 5 (Data Representation Section) of a GRIB2 message, starting at octet 6 of that Sectio...
Definition: gf_unpack5.F90:35
Handles Data Representation Templates used in Section 5.
subroutine getdrstemplate(number, nummap, map, needext, iret)
Return DRS template information for a specified Data Representation Template.
subroutine extdrstemplate(number, list, nummap, map)
Generate the remaining octet map for a given Data Representation Template, if required.
Reallocate memory, preserving contents.
Definition: realloc.f:12