NCEPLIBS-g2  3.4.7
compack.f
Go to the documentation of this file.
1 
4 
34  subroutine compack(fld,ndpts,idrsnum,idrstmpl,cpack,lcpack)
35 
36  use intmath
37  implicit none
38  integer,intent(in) :: ndpts,idrsnum
39  real,intent(in) :: fld(ndpts)
40  character(len=1),intent(out) :: cpack(*)
41  integer,intent(inout) :: idrstmpl(*)
42  integer,intent(out) :: lcpack
43 
44  real :: bscale,dscale
45  integer :: j,iofst,imin,ival1,ival2,minsd,nbitsd,n
46  integer :: igmax,nbitsgref,left,iwmax,i,ilmax,kk,ij
47  integer :: ngwidthref,nbitsgwidth,nglenref,nglenlast
48  integer :: maxorig,nbitorig,isd,ngroups,itemp,minpk
49  integer :: kfildo,inc,maxgrps,missopt,miss1,miss2,lg
50  integer :: ibit,jbit,kbit,novref,lbitref,ier,ng,imax
51  integer :: nbitsglen
52  real(4) :: ref,rmin4
53  real(8) :: rmin,rmax
54 
55  integer(4) :: iref
56  integer,allocatable :: ifld(:)
57  integer,allocatable :: jmin(:),jmax(:),lbit(:)
58  integer,parameter :: zero=0
59  integer,allocatable :: gref(:),gwidth(:),glen(:)
60  integer :: glength,grpwidth
61  logical :: simple_alg
62 
63  simple_alg = .false.
64 
65  bscale=2.0**real(-idrstmpl(2))
66  dscale=10.0**real(idrstmpl(3))
67 !
68 ! Find max and min values in the data
69 !
70  if(ndpts>0) then
71  rmax=fld(1)
72  rmin=fld(1)
73  else
74  rmax=1.0
75  rmin=1.0
76  endif
77  do j=2,ndpts
78  if (fld(j).gt.rmax) rmax=fld(j)
79  if (fld(j).lt.rmin) rmin=fld(j)
80  enddo
81 !
82 ! If max and min values are not equal, pack up field.
83 ! If they are equal, we have a constant field, and the reference
84 ! value (rmin) is the value for each point in the field and
85 ! set nbits to 0.
86 !
87  multival: if (rmin.ne.rmax) then
88  iofst=0
89  allocate(ifld(ndpts))
90  allocate(gref(ndpts))
91  allocate(gwidth(ndpts))
92  allocate(glen(ndpts))
93  !
94  ! Scale original data
95  !
96  if (idrstmpl(2).eq.0) then ! No binary scaling
97  imin=nint(rmin*dscale)
98  !imax=nint(rmax*dscale)
99  rmin=real(imin)
100  do j=1,ndpts
101  ifld(j)=max(0,nint(fld(j)*dscale)-imin)
102  enddo
103  else ! Use binary scaling factor
104  rmin=rmin*dscale
105  !rmax=rmax*dscale
106  do j=1,ndpts
107  ifld(j)=max(0,nint(((fld(j)*dscale)-rmin)*bscale))
108  enddo
109  endif
110  !
111  ! Calculate Spatial differences, if using DRS Template 5.3
112  !
113  alg3: if (idrsnum.eq.3) then ! spatial differences
114  if (idrstmpl(17).ne.1.and.idrstmpl(17).ne.2) idrstmpl(17)=2
115  if (idrstmpl(17).eq.1) then ! first order
116  ival1=ifld(1)
117  if(ival1<0) then
118  print *,'G2: negative ival1',ival1
119  stop 101
120  endif
121  do j=ndpts,2,-1
122  ifld(j)=ifld(j)-ifld(j-1)
123  enddo
124  ifld(1)=0
125  elseif (idrstmpl(17).eq.2) then ! second order
126  ival1=ifld(1)
127  ival2=ifld(2)
128  if(ival1<0) then
129  print *,'G2: negative ival1',ival1
130  stop 11
131  endif
132  if(ival2<0) then
133  print *,'G2: negative ival2',ival2
134  stop 12
135  endif
136  do j=ndpts,3,-1
137  ifld(j)=ifld(j)-(2*ifld(j-1))+ifld(j-2)
138  enddo
139  ifld(1)=0
140  ifld(2)=0
141  endif
142  !
143  ! subtract min value from spatial diff field
144  !
145  isd=idrstmpl(17)+1
146  minsd=minval(ifld(isd:ndpts))
147  do j=isd,ndpts
148  ifld(j)=ifld(j)-minsd
149  enddo
150  !
151  ! find num of bits need to store minsd and add 1 extra bit
152  ! to indicate sign
153  !
154  nbitsd=i1log2(abs(minsd))+1
155  !
156  ! find num of bits need to store ifld(1) ( and ifld(2)
157  ! if using 2nd order differencing )
158  !
159  maxorig=ival1
160  if (idrstmpl(17).eq.2.and.ival2.gt.ival1) maxorig=ival2
161  nbitorig=i1log2(maxorig)+1
162  if (nbitorig.gt.nbitsd) nbitsd=nbitorig
163  ! increase number of bits to even multiple of 8 ( octet )
164  if (mod(nbitsd,8).ne.0) nbitsd=nbitsd+(8-mod(nbitsd,8))
165  !
166  ! Store extra spatial differencing info into the packed
167  ! data section.
168  !
169  if (nbitsd.ne.0) then
170  ! pack first original value
171  if (ival1.ge.0) then
172  call g2_sbytec(cpack,ival1,iofst,nbitsd)
173  iofst=iofst+nbitsd
174  else
175  call g2_sbytec(cpack,1,iofst,1)
176  iofst=iofst+1
177  call g2_sbytec(cpack,iabs(ival1),iofst,nbitsd-1)
178  iofst=iofst+nbitsd-1
179  endif
180  if (idrstmpl(17).eq.2) then
181  ! pack second original value
182  if (ival2.ge.0) then
183  call g2_sbytec(cpack,ival2,iofst,nbitsd)
184  iofst=iofst+nbitsd
185  else
186  call g2_sbytec(cpack,1,iofst,1)
187  iofst=iofst+1
188  call g2_sbytec(cpack,iabs(ival2),iofst,nbitsd-1)
189  iofst=iofst+nbitsd-1
190  endif
191  endif
192  ! pack overall min of spatial differences
193  if (minsd.ge.0) then
194  call g2_sbytec(cpack,minsd,iofst,nbitsd)
195  iofst=iofst+nbitsd
196  else
197  call g2_sbytec(cpack,1,iofst,1)
198  iofst=iofst+1
199  call g2_sbytec(cpack,iabs(minsd),iofst,nbitsd-1)
200  iofst=iofst+nbitsd-1
201  endif
202  endif
203  !print *,'SDp ',ival1,ival2,minsd,nbitsd
204  endif alg3 ! end of spatial diff section
205  !
206  ! Determine Groups to be used.
207  !
208  simplealg: if ( simple_alg ) then
209  ! set group length to 10 : calculate number of groups
210  ! and length of last group
211  print *,'G2: use simple algorithm'
212  ngroups=ndpts/10
213  glen(1:ngroups)=10
214  itemp=mod(ndpts,10)
215  if (itemp.ne.0) then
216  ngroups=ngroups+1
217  glen(ngroups)=itemp
218  endif
219  else
220  ! Use Dr. Glahn's algorithm for determining grouping.
221  !
222  kfildo=6
223  minpk=10
224  inc=1
225  maxgrps=((ndpts+minpk-1)/minpk)
226  allocate(jmin(maxgrps))
227  allocate(jmax(maxgrps))
228  allocate(lbit(maxgrps))
229  missopt=0
230  call pack_gp(kfildo,ifld,ndpts,missopt,minpk,inc,miss1,miss2,
231  & jmin,jmax,lbit,glen,maxgrps,ngroups,ibit,jbit,
232  & kbit,novref,lbitref,ier)
233  if(ier/=0) then
234  ! Dr. Glahn's algorithm failed; use simple packing method instead.
235  1099 format('G2: fall back to simple algorithm (glahn ier=',i0,&
236  & ')')
237  print 1099,ier
238  ngroups=ndpts/10
239  glen(1:ngroups)=10
240  itemp=mod(ndpts,10)
241  if (itemp.ne.0) then
242  ngroups=ngroups+1
243  glen(ngroups)=itemp
244  endif
245  elseif(ngroups<1) then
246  ! Dr. Glahn's algorithm failed; use simple packing method instead.
247  print *,'Glahn algorithm failed; use simple packing'
248  ngroups=ndpts/10
249  glen(1:ngroups)=10
250  itemp=mod(ndpts,10)
251  if (itemp.ne.0) then
252  ngroups=ngroups+1
253  glen(ngroups)=itemp
254  endif
255  else
256 !print *,'SAGier = ',ier,ibit,jbit,kbit,novref,lbitref
257  do ng=1,ngroups
258  glen(ng)=glen(ng)+novref
259  enddo
260  deallocate(jmin)
261  deallocate(jmax)
262  deallocate(lbit)
263  endif
264  endif simplealg
265  !
266  ! For each group, find the group's reference value
267  ! and the number of bits needed to hold the remaining values
268  !
269  n=1
270  do ng=1,ngroups
271  ! find max and min values of group
272  gref(ng)=ifld(n)
273  imax=ifld(n)
274  j=n+1
275  do lg=2,glen(ng)
276  if (ifld(j).lt.gref(ng)) gref(ng)=ifld(j)
277  if (ifld(j).gt.imax) imax=ifld(j)
278  j=j+1
279  enddo
280  ! calc num of bits needed to hold data
281  if ( gref(ng).ne.imax ) then
282  gwidth(ng)=i1log2(imax-gref(ng))
283  else
284  gwidth(ng)=0
285  endif
286  ! Subtract min from data
287  j=n
288  do lg=1,glen(ng)
289  ifld(j)=ifld(j)-gref(ng)
290  j=j+1
291  enddo
292  ! increment fld array counter
293  n=n+glen(ng)
294  enddo
295  !
296  ! Find max of the group references and calc num of bits needed
297  ! to pack each groups reference value, then
298  ! pack up group reference values
299  !
300  !write(77,*)'GREFS: ',(gref(j),j=1,ngroups)
301  igmax=maxval(gref(1:ngroups))
302  if (igmax.ne.0) then
303  nbitsgref=i1log2(igmax)
304  call g2_sbytesc(cpack,gref,iofst,nbitsgref,0,ngroups)
305  itemp=nbitsgref*ngroups
306  iofst=iofst+itemp
307  if (mod(itemp,8).ne.0) then
308  left=8-mod(itemp,8)
309  call g2_sbytec(cpack,zero,iofst,left)
310  iofst=iofst+left
311  endif
312  else
313  nbitsgref=0
314  endif
315  !
316  ! Find max/min of the group widths and calc num of bits needed
317  ! to pack each groups width value, then
318  ! pack up group width values
319  !
320  !write(77,*)'GWIDTHS: ',(gwidth(j),j=1,ngroups)
321  iwmax=maxval(gwidth(1:ngroups))
322  ngwidthref=minval(gwidth(1:ngroups))
323  if (iwmax.ne.ngwidthref) then
324  nbitsgwidth=i1log2(iwmax-ngwidthref)
325  do i=1,ngroups
326  gwidth(i)=gwidth(i)-ngwidthref
327  if(gwidth(i)<0) then
328  write(0,*) 'i,gw,ngw=',i,gwidth(i),ngwidthref
329  stop 9
330  endif
331  enddo
332  call g2_sbytesc(cpack,gwidth,iofst,nbitsgwidth,0,ngroups)
333  itemp=nbitsgwidth*ngroups
334  iofst=iofst+itemp
335  ! Pad last octet with Zeros, if necessary,
336  if (mod(itemp,8).ne.0) then
337  left=8-mod(itemp,8)
338  call g2_sbytec(cpack,zero,iofst,left)
339  iofst=iofst+left
340  endif
341  else
342  nbitsgwidth=0
343  gwidth(1:ngroups)=0
344  endif
345  !
346  ! Find max/min of the group lengths and calc num of bits needed
347  ! to pack each groups length value, then
348  ! pack up group length values
349  !
350  !write(77,*)'GLENS: ',(glen(j),j=1,ngroups)
351  ilmax=maxval(glen(1:ngroups-1))
352  nglenref=minval(glen(1:ngroups-1))
353  nglenlast=glen(ngroups)
354  if (ilmax.ne.nglenref) then
355  nbitsglen=i1log2(ilmax-nglenref)
356  do i=1,ngroups-1
357  glen(i)=glen(i)-nglenref
358  if(glen(i)<0) then
359  write(0,*) 'i,glen(i) = ',i,glen(i)
360  stop 23
361  endif
362  enddo
363  call g2_sbytesc(cpack,glen,iofst,nbitsglen,0,ngroups)
364  itemp=nbitsglen*ngroups
365  iofst=iofst+itemp
366  ! Pad last octet with Zeros, if necessary,
367  if (mod(itemp,8).ne.0) then
368  left=8-mod(itemp,8)
369  call g2_sbytec(cpack,zero,iofst,left)
370  iofst=iofst+left
371  endif
372  else
373  nbitsglen=0
374  glen(1:ngroups)=0
375  endif
376  !
377  ! For each group, pack data values
378  !
379  !write(77,*)'IFLDS: ',(ifld(j),j=1,ndpts)
380  n=1
381  ij=0
382  do ng=1,ngroups
383  glength=glen(ng)+nglenref
384  if (ng.eq.ngroups ) glength=nglenlast
385  grpwidth=gwidth(ng)+ngwidthref
386  !write(77,*)'NGP ',ng,grpwidth,glength,gref(ng)
387  if ( grpwidth.ne.0 ) then
388  call g2_sbytesc(cpack,ifld(n),iofst,grpwidth,0,glength)
389  iofst=iofst+(grpwidth*glength)
390  endif
391  do kk=1,glength
392  ij=ij+1
393  !write(77,*)'SAG ',ij,fld(ij),ifld(ij),gref(ng),bscale,rmin,dscale
394  enddo
395  n=n+glength
396  enddo
397  ! Pad last octet with Zeros, if necessary,
398  if (mod(iofst,8).ne.0) then
399  left=8-mod(iofst,8)
400  call g2_sbytec(cpack,zero,iofst,left)
401  iofst=iofst+left
402  endif
403  lcpack=iofst/8
404  !
405  if ( allocated(ifld) ) deallocate(ifld)
406  if ( allocated(gref) ) deallocate(gref)
407  if ( allocated(gwidth) ) deallocate(gwidth)
408  if ( allocated(glen) ) deallocate(glen)
409  else ! Constant field ( max = min )
410  lcpack=0
411  nbitsgref=0
412  ngroups=0
413  ngwidthref=0
414  nbitsgwidth=0
415  nglenref=0
416  nglenlast=0
417  nbitsglen=0
418  nbitsd=0
419  endif multival
420 
421 !
422 ! Fill in ref value and number of bits in Template 5.2
423 !
424  rmin4 = real(rmin, 4)
425  call mkieee(rmin4,ref,1) ! ensure reference value is IEEE format
426 ! call g2_gbytec(ref,idrstmpl(1),0,32)
427  iref=transfer(ref,iref)
428  idrstmpl(1)=iref
429  idrstmpl(4)=nbitsgref
430  idrstmpl(5)=0 ! original data were reals
431  idrstmpl(6)=1 ! general group splitting
432  idrstmpl(7)=0 ! No internal missing values
433  idrstmpl(8)=0 ! Primary missing value
434  idrstmpl(9)=0 ! secondary missing value
435  idrstmpl(10)=ngroups ! Number of groups
436  idrstmpl(11)=ngwidthref ! reference for group widths
437  idrstmpl(12)=nbitsgwidth ! num bits used for group widths
438  idrstmpl(13)=nglenref ! Reference for group lengths
439  idrstmpl(14)=1 ! length increment for group lengths
440  idrstmpl(15)=nglenlast ! True length of last group
441  idrstmpl(16)=nbitsglen ! num bits used for group lengths
442  if (idrsnum.eq.3) then
443  idrstmpl(18)=nbitsd/8 ! num bits used for extra spatial
444  ! differencing values
445  endif
446 
447  end
subroutine compack(fld, ndpts, idrsnum, idrstmpl, cpack, lcpack)
This subroutine supports GRIB2 complex packing templates with or without spatial differences,...
Definition: compack.f:35
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 g2_sbytesc(out, in, iskip, nbits, nskip, n)
Put arbitrary size values into a packed bit string, taking the low order bits from each value in the ...
Definition: g2_gbytesc.F90:120
subroutine mkieee(a, rieee, num)
Copy an array of real to an array of 32-bit IEEE floating points.
Definition: mkieee.F90:14
Define math functions used by compack(), simpack(), and misspack().
Definition: intmath.f:22
subroutine pack_gp(KFILDO, IC, NXY, IS523, MINPK, INC, MISSP, MISSS, JMIN, JMAX, LBIT, NOV, NDG, LX, IBIT, JBIT, KBIT, NOVREF, LBITREF, IER)
This subroutine determines groups of variable size, but at least of size minpk, the associated max(JM...
Definition: pack_gp.f:117