NCEPLIBS-g2  3.4.5
gf_unpack5.f
Go to the documentation of this file.
1 
6 
35 
36  subroutine gf_unpack5(cgrib,lcgrib,iofst,ndpts,idrsnum,idrstmpl,
37  & mapdrslen,ierr)
38 
39  use drstemplates
40  use re_alloc ! needed for subroutine realloc
41 
42  character(len=1),intent(in) :: cgrib(lcgrib)
43  integer,intent(in) :: lcgrib
44  integer,intent(inout) :: iofst
45  integer,intent(out) :: ndpts,idrsnum
46  integer,pointer,dimension(:) :: idrstmpl
47  integer,intent(out) :: ierr
48 
49  integer,allocatable :: mapdrs(:)
50  integer :: mapdrslen
51  logical needext
52 
53  ierr=0
54  nullify(idrstmpl)
55 
56  call g2_gbytec(cgrib,lensec,iofst,32) ! Get Length of Section
57  iofst=iofst+32
58  iofst=iofst+8 ! skip section number
59  allocate(mapdrs(lensec))
60 
61  call g2_gbytec(cgrib,ndpts,iofst,32) ! Get num of data points
62  iofst=iofst+32
63  call g2_gbytec(cgrib,idrsnum,iofst,16) ! Get Data Rep Template Num.
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
74  ! the appropriate number of octets, which are specified in
75  ! corresponding entries in array mappds.
76  !
77  istat=0
78  if (mapdrslen.gt.0) allocate(idrstmpl(mapdrslen),stat=istat)
79  if (istat.ne.0) then
80  ierr=6
81  nullify(idrstmpl)
82  if( allocated(mapdrs) ) deallocate(mapdrs)
83  return
84  endif
85  do i=1,mapdrslen
86  nbits=iabs(mapdrs(i))*8
87  if ( mapdrs(i).ge.0 ) then
88  call g2_gbytec(cgrib,idrstmpl(i),iofst,nbits)
89  else
90  call g2_gbytec(cgrib,isign,iofst,1)
91  call g2_gbytec(cgrib,idrstmpl(i),iofst+1,nbits-1)
92  if (isign.eq.1) idrstmpl(i)=-idrstmpl(i)
93  endif
94  iofst=iofst+nbits
95  enddo
96  !
97  ! Check to see if the Data Representation Template needs to be
98  ! extended.
99  ! The number of values in a specific template may vary
100  ! depending on data specified in the "static" part of the
101  ! template.
102  !
103  if ( needext ) then
104  call extdrstemplate(idrsnum,idrstmpl,newmapdrslen,mapdrs)
105  call realloc(idrstmpl,mapdrslen,newmapdrslen,istat)
106  ! Unpack the rest of the Data Representation Template
107  do i=mapdrslen+1,newmapdrslen
108  nbits=iabs(mapdrs(i))*8
109  if ( mapdrs(i).ge.0 ) then
110  call g2_gbytec(cgrib,idrstmpl(i),iofst,nbits)
111  else
112  call g2_gbytec(cgrib,isign,iofst,1)
113  call g2_gbytec(cgrib,idrstmpl(i),iofst+1,nbits-1)
114  if (isign.eq.1) idrstmpl(i)=-idrstmpl(i)
115  endif
116  iofst=iofst+nbits
117  enddo
118  mapdrslen=newmapdrslen
119  endif
120  if( allocated(mapdrs) ) deallocate(mapdrs)
121 
122  return ! End of Section 5 processing
123  end
124 
drstemplates
This Fortran Module contains info on all the available GRIB2 Data Representation Templates used in Se...
Definition: drstemplates.f:37
drstemplates::getdrstemplate
subroutine getdrstemplate(number, nummap, map, needext, iret)
This subroutine returns DRS template information for a .
Definition: drstemplates.f:165
re_alloc
This module contains three subroutines to reorganize the integer, real and character data in memory i...
Definition: realloc.f:14
re_alloc::realloc
Definition: realloc.f:16
drstemplates::extdrstemplate
subroutine extdrstemplate(number, list, nummap, map)
This subroutine generates the remaining octet map for a given Data Representation Template,...
Definition: drstemplates.f:207
gf_unpack5
subroutine gf_unpack5(cgrib, lcgrib, iofst, ndpts, idrsnum, idrstmpl, mapdrslen, ierr)
This subroutine unpacks Section 5 (Data Representation Section) starting at octet 6 of that Section.
Definition: gf_unpack5.f:38
g2_gbytec
subroutine g2_gbytec(IN, IOUT, ISKIP, NBYTE)
This subrountine is to extract arbitrary size values from a packed bit string, right justifying each ...
Definition: g2_gbytesc.f:20