NCEPLIBS-g2  3.4.5
misspack.f
Go to the documentation of this file.
1 
5 
42 
43  subroutine misspack(fld,ndpts,idrsnum,idrstmpl,cpack,lcpack)
44 
45  use intmath
46  integer,intent(in) :: ndpts,idrsnum
47  real,intent(in) :: fld(ndpts)
48  character(len=1),intent(out) :: cpack(*)
49  integer,intent(inout) :: idrstmpl(*)
50  integer,intent(out) :: lcpack
51 
52  real(4) :: ref, rmin4
53  integer(4) :: iref
54  integer,allocatable :: ifld(:),ifldmiss(:),jfld(:)
55  integer,allocatable :: jmin(:),jmax(:),lbit(:)
56  integer,parameter :: zero=0
57  integer,allocatable :: gref(:),gwidth(:),glen(:)
58  integer :: glength,grpwidth
59  logical :: simple_alg
60  logical :: have_rmin
61 
62  simple_alg = .false.
63  have_rmin = .false.
64  bscale=2.0**real(-idrstmpl(2))
65  dscale=10.0**real(idrstmpl(3))
66  missopt=idrstmpl(7)
67  if ( missopt.ne.1 .AND. missopt.ne.2 ) then
68  print *,'misspack: Unrecognized option.'
69  lcpack=-1
70  return
71  else ! Get missing values
72  call rdieee(idrstmpl(8),rmissp,1)
73  if (missopt.eq.2) call rdieee(idrstmpl(9),rmisss,1)
74  endif
75 !
76 ! Find min value of non-missing values in the data,
77 ! AND set up missing value mapping of the field.
78 !
79  allocate(ifldmiss(ndpts))
80 c rmin=huge(rmin)
81 
82  if ( missopt .eq. 1 ) then ! Primary missing value only
83  do j=1,ndpts
84  if (fld(j).eq.rmissp) then
85  ifldmiss(j)=1
86  else
87  ifldmiss(j)=0
88  if(have_rmin) then
89  if (fld(j).lt.rmin) rmin=fld(j)
90  else
91  rmin=fld(j)
92  have_rmin=.true.
93  endif
94  endif
95  enddo
96  if(.not.have_rmin) rmin=rmissp
97  endif
98  if ( missopt .eq. 2 ) then ! Primary and secondary missing values
99  do j=1,ndpts
100  if (fld(j).eq.rmissp) then
101  ifldmiss(j)=1
102  elseif (fld(j).eq.rmisss) then
103  ifldmiss(j)=2
104  else
105  ifldmiss(j)=0
106  if(have_rmin) then
107  if (fld(j).lt.rmin) rmin=fld(j)
108  else
109  rmin=fld(j)
110  have_rmin=.true.
111  endif
112  endif
113  if(.not.have_rmin) rmin=rmissp
114  enddo
115  endif
116 !
117 ! Allocate work arrays:
118 ! Note: -ifldmiss(j),j=1,ndpts is a map of original field indicating
119 ! which of the original data values
120 ! are primary missing (1), sencondary missing (2) or non-missing (0).
121 ! -jfld(j),j=1,nonmiss is a subarray of just the non-missing values from
122 ! the original field.
123 !
124  !if (rmin.ne.rmax) then
125  iofst=0
126  allocate(ifld(ndpts))
127  allocate(jfld(ndpts))
128  allocate(gref(ndpts))
129  allocate(gwidth(ndpts))
130  allocate(glen(ndpts))
131  !
132  ! Scale original data
133  !
134  nonmiss=0
135  if (idrstmpl(2).eq.0) then ! No binary scaling
136  imin=nint(rmin*dscale)
137  !imax=nint(rmax*dscale)
138  rmin=real(imin)
139  do j=1,ndpts
140  if (ifldmiss(j).eq.0) then
141  nonmiss=nonmiss+1
142  jfld(nonmiss)=max(0,nint(fld(j)*dscale)-imin)
143  endif
144  enddo
145  else ! Use binary scaling factor
146  rmin=rmin*dscale
147  !rmax=rmax*dscale
148  do j=1,ndpts
149  if (ifldmiss(j).eq.0) then
150  nonmiss=nonmiss+1
151  jfld(nonmiss)=max(0,nint(((fld(j)*dscale)-rmin)*bscale))
152  endif
153  enddo
154  endif
155  !
156  ! Calculate Spatial differences, if using DRS Template 5.3
157  !
158  if (idrsnum.eq.3) then ! spatial differences
159  if (idrstmpl(17).ne.1.and.idrstmpl(17).ne.2) idrstmpl(17)=2
160  if (idrstmpl(17).eq.1) then ! first order
161  if(nonmiss<1) then
162  ival1=1.0
163  else
164  ival1=jfld(1)
165  endif
166  do j=nonmiss,2,-1
167  jfld(j)=jfld(j)-jfld(j-1)
168  enddo
169  if(nonmiss>0) jfld(1)=0
170  elseif (idrstmpl(17).eq.2) then ! second order
171  if(nonmiss==1) then
172  ival1=jfld(1)
173  ival2=jfld(1)
174  elseif(nonmiss<1) then
175  ival1=1.0
176  ival2=1.0
177  else
178  ival1=jfld(1)
179  ival2=jfld(2)
180  endif
181  do j=nonmiss,3,-1
182  jfld(j)=jfld(j)-(2*jfld(j-1))+jfld(j-2)
183  enddo
184  if(nonmiss>=1) jfld(1)=0
185  if(nonmiss>=2) jfld(2)=0
186  endif
187  !
188  ! subtract min value from spatial diff field
189  !
190  isd=idrstmpl(17)+1
191  minsd=minval(jfld(isd:nonmiss))
192  do j=isd,nonmiss
193  jfld(j)=jfld(j)-minsd
194  enddo
195  !
196  ! find num of bits need to store minsd and add 1 extra bit
197  ! to indicate sign
198  !
199  temp=i1log2(abs(minsd))
200  nbitsd=ceiling(temp)+1
201  !
202  ! find num of bits need to store ifld(1) ( and ifld(2)
203  ! if using 2nd order differencing )
204  !
205  maxorig=ival1
206  if (idrstmpl(17).eq.2.and.ival2.gt.ival1) maxorig=ival2
207  temp=i1log2(maxorig)
208  nbitorig=ceiling(temp)+1
209  if (nbitorig.gt.nbitsd) nbitsd=nbitorig
210  ! increase number of bits to even multiple of 8 ( octet )
211  if (mod(nbitsd,8).ne.0) nbitsd=nbitsd+(8-mod(nbitsd,8))
212  !
213  ! Store extra spatial differencing info into the packed
214  ! data section.
215  !
216  if (nbitsd.ne.0) then
217  ! pack first original value
218  if (ival1.ge.0) then
219  call g2_sbytec(cpack,ival1,iofst,nbitsd)
220  iofst=iofst+nbitsd
221  else
222  call g2_sbytec(cpack,1,iofst,1)
223  iofst=iofst+1
224  call g2_sbytec(cpack,iabs(ival1),iofst,nbitsd-1)
225  iofst=iofst+nbitsd-1
226  endif
227  if (idrstmpl(17).eq.2) then
228  ! pack second original value
229  if (ival2.ge.0) then
230  call g2_sbytec(cpack,ival2,iofst,nbitsd)
231  iofst=iofst+nbitsd
232  else
233  call g2_sbytec(cpack,1,iofst,1)
234  iofst=iofst+1
235  call g2_sbytec(cpack,iabs(ival2),iofst,nbitsd-1)
236  iofst=iofst+nbitsd-1
237  endif
238  endif
239  ! pack overall min of spatial differences
240  if (minsd.ge.0) then
241  call g2_sbytec(cpack,minsd,iofst,nbitsd)
242  iofst=iofst+nbitsd
243  else
244  call g2_sbytec(cpack,1,iofst,1)
245  iofst=iofst+1
246  call g2_sbytec(cpack,iabs(minsd),iofst,nbitsd-1)
247  iofst=iofst+nbitsd-1
248  endif
249  endif
250  !print *,'SDp ',ival1,ival2,minsd,nbitsd
251  endif ! end of spatial diff section
252  !
253  ! Expand non-missing data values to original grid.
254  !
255  miss1=minval(jfld(1:nonmiss))-1
256  miss2=miss1-1
257  n=0
258  do j=1,ndpts
259  if ( ifldmiss(j).eq.0 ) then
260  n=n+1
261  ifld(j)=jfld(n)
262  elseif ( ifldmiss(j).eq.1 ) then
263  ifld(j)=miss1
264  elseif ( ifldmiss(j).eq.2 ) then
265  ifld(j)=miss2
266  endif
267  enddo
268  if(ndpts<2) simple_alg=.true.
269  !
270  ! Determine Groups to be used.
271  !
272  if ( simple_alg ) then
273  ! set group length to 10 : calculate number of groups
274  ! and length of last group
275  ngroups=ndpts/10
276  glen(1:ngroups)=10
277  itemp=mod(ndpts,10)
278  if (itemp.ne.0) then
279  ngroups=ngroups+1
280  glen(ngroups)=itemp
281  endif
282  else
283  ! Use Dr. Glahn's algorithm for determining grouping.
284  !
285  kfildo=6
286  minpk=10
287  inc=1
288  maxgrps=(ndpts/minpk)+1
289  allocate(jmin(maxgrps))
290  allocate(jmax(maxgrps))
291  allocate(lbit(maxgrps))
292  call pack_gp(kfildo,ifld,ndpts,missopt,minpk,inc,miss1,miss2,
293  & jmin,jmax,lbit,glen,maxgrps,ngroups,ibit,jbit,
294  & kbit,novref,lbitref,ier)
295  !print *,'SAGier = ',ier,ibit,jbit,kbit,novref,lbitref
296  do ng=1,ngroups
297  glen(ng)=glen(ng)+novref
298  enddo
299  deallocate(jmin)
300  deallocate(jmax)
301  deallocate(lbit)
302  endif
303  !
304  ! For each group, find the group's reference value (min)
305  ! and the number of bits needed to hold the remaining values
306  !
307  n=1
308  do ng=1,ngroups
309  ! how many of each type?
310  num0=count(ifldmiss(n:n+glen(ng)-1) .EQ. 0)
311  num1=count(ifldmiss(n:n+glen(ng)-1) .EQ. 1)
312  num2=count(ifldmiss(n:n+glen(ng)-1) .EQ. 2)
313  if ( num0.eq.0 ) then ! all missing values
314  if ( num1.eq.0 ) then ! all secondary missing
315  gref(ng)=-2
316  gwidth(ng)=0
317  elseif ( num2.eq.0 ) then ! all primary missing
318  gref(ng)=-1
319  gwidth(ng)=0
320  else ! both primary and secondary
321  gref(ng)=0
322  gwidth(ng)=1
323  endif
324  else ! contains some non-missing data
325  ! find max and min values of group
326  gref(ng)=huge(n)
327  imax=-1*huge(n)
328  j=n
329  do lg=1,glen(ng)
330  if ( ifldmiss(j).eq.0 ) then
331  if (ifld(j).lt.gref(ng)) gref(ng)=ifld(j)
332  if (ifld(j).gt.imax) imax=ifld(j)
333  endif
334  j=j+1
335  enddo
336  if (missopt.eq.1) imax=imax+1
337  if (missopt.eq.2) imax=imax+2
338  ! calc num of bits needed to hold data
339  if ( gref(ng).ne.imax ) then
340  temp=i1log2(imax-gref(ng))
341  gwidth(ng)=ceiling(temp)
342  else
343  gwidth(ng)=0
344  endif
345  endif
346  ! Subtract min from data
347  j=n
348  mtemp=2**gwidth(ng)
349  do lg=1,glen(ng)
350  if (ifldmiss(j).eq.0) then ! non-missing
351  ifld(j)=ifld(j)-gref(ng)
352  elseif (ifldmiss(j).eq.1) then ! primary missing
353  ifld(j)=mtemp-1
354  elseif (ifldmiss(j).eq.2) then ! secondary missing
355  ifld(j)=mtemp-2
356  endif
357  j=j+1
358  enddo
359  ! increment fld array counter
360  n=n+glen(ng)
361  enddo
362  !
363  ! Find max of the group references and calc num of bits needed
364  ! to pack each groups reference value, then
365  ! pack up group reference values
366  !
367  !write(77,*)'GREFS: ',(gref(j),j=1,ngroups)
368  igmax=maxval(gref(1:ngroups))
369  if (missopt.eq.1) igmax=igmax+1
370  if (missopt.eq.2) igmax=igmax+2
371  if (igmax.ne.0) then
372  temp=i1log2(igmax)
373  nbitsgref=ceiling(temp)
374  ! restet the ref values of any "missing only" groups.
375  mtemp=2**nbitsgref
376  do j=1,ngroups
377  if (gref(j).eq.-1) gref(j)=mtemp-1
378  if (gref(j).eq.-2) gref(j)=mtemp-2
379  enddo
380  call g2_sbytesc(cpack,gref,iofst,nbitsgref,0,ngroups)
381  itemp=nbitsgref*ngroups
382  iofst=iofst+itemp
383  ! Pad last octet with Zeros, if necessary,
384  if (mod(itemp,8).ne.0) then
385  left=8-mod(itemp,8)
386  call g2_sbytec(cpack,zero,iofst,left)
387  iofst=iofst+left
388  endif
389  else
390  nbitsgref=0
391  endif
392  !
393  ! Find max/min of the group widths and calc num of bits needed
394  ! to pack each groups width value, then
395  ! pack up group width values
396  !
397  !write(77,*)'GWIDTHS: ',(gwidth(j),j=1,ngroups)
398  iwmax=maxval(gwidth(1:ngroups))
399  ngwidthref=minval(gwidth(1:ngroups))
400  if (iwmax.ne.ngwidthref) then
401  temp=i1log2(iwmax-ngwidthref)
402  nbitsgwidth=ceiling(temp)
403  do i=1,ngroups
404  gwidth(i)=gwidth(i)-ngwidthref
405  enddo
406  call g2_sbytesc(cpack,gwidth,iofst,nbitsgwidth,0,ngroups)
407  itemp=nbitsgwidth*ngroups
408  iofst=iofst+itemp
409  ! Pad last octet with Zeros, if necessary,
410  if (mod(itemp,8).ne.0) then
411  left=8-mod(itemp,8)
412  call g2_sbytec(cpack,zero,iofst,left)
413  iofst=iofst+left
414  endif
415  else
416  nbitsgwidth=0
417  gwidth(1:ngroups)=0
418  endif
419  !
420  ! Find max/min of the group lengths and calc num of bits needed
421  ! to pack each groups length value, then
422  ! pack up group length values
423  !
424  !write(77,*)'GLENS: ',(glen(j),j=1,ngroups)
425  ilmax=maxval(glen(1:ngroups-1))
426  nglenref=minval(glen(1:ngroups-1))
427  if(ngroups>0) then
428  nglenlast=glen(ngroups)
429  else
430  nglenlast=0
431  endif
432  if (ilmax.ne.nglenref) then
433  temp=i1log2(ilmax-nglenref)
434  nbitsglen=ceiling(temp)
435  do i=1,ngroups-1
436  glen(i)=glen(i)-nglenref
437  enddo
438  call g2_sbytesc(cpack,glen,iofst,nbitsglen,0,ngroups)
439  itemp=nbitsglen*ngroups
440  iofst=iofst+itemp
441  ! Pad last octet with Zeros, if necessary,
442  if (mod(itemp,8).ne.0) then
443  left=8-mod(itemp,8)
444  call g2_sbytec(cpack,zero,iofst,left)
445  iofst=iofst+left
446  endif
447  else
448  nbitsglen=0
449  glen(1:ngroups)=0
450  endif
451  !
452  ! For each group, pack data values
453  !
454  !write(77,*)'IFLDS: ',(ifld(j),j=1,ndpts)
455  n=1
456  ij=0
457  do ng=1,ngroups
458  glength=glen(ng)+nglenref
459  if (ng.eq.ngroups ) glength=nglenlast
460  grpwidth=gwidth(ng)+ngwidthref
461  !write(77,*)'NGP ',ng,grpwidth,glength,gref(ng)
462  if ( grpwidth.ne.0 ) then
463  call g2_sbytesc(cpack,ifld(n),iofst,grpwidth,0,glength)
464  iofst=iofst+(grpwidth*glength)
465  endif
466  do kk=1,glength
467  ij=ij+1
468  !write(77,*)'SAG ',ij,fld(ij),ifld(ij),gref(ng),bscale,rmin,dscale
469  enddo
470  n=n+glength
471  enddo
472  ! Pad last octet with Zeros, if necessary,
473  if (mod(iofst,8).ne.0) then
474  left=8-mod(iofst,8)
475  call g2_sbytec(cpack,zero,iofst,left)
476  iofst=iofst+left
477  endif
478  lcpack=iofst/8
479  !
480  if ( allocated(ifld) ) deallocate(ifld)
481  if ( allocated(jfld) ) deallocate(jfld)
482  if ( allocated(ifldmiss) ) deallocate(ifldmiss)
483  if ( allocated(gref) ) deallocate(gref)
484  if ( allocated(gwidth) ) deallocate(gwidth)
485  if ( allocated(glen) ) deallocate(glen)
486  !else ! Constant field ( max = min )
487  ! nbits=0
488  ! lcpack=0
489  ! nbitsgref=0
490  ! ngroups=0
491  !endif
492 
493 !
494 ! Fill in ref value and number of bits in Template 5.2
495 !
496  rmin4 = rmin
497  call mkieee(rmin4,ref,1) ! ensure reference value is IEEE format
498 ! call g2_gbytec(ref,idrstmpl(1),0,32)
499  iref=transfer(ref,iref)
500  idrstmpl(1)=iref
501  idrstmpl(4)=nbitsgref
502  idrstmpl(5)=0 ! original data were reals
503  idrstmpl(6)=1 ! general group splitting
504  idrstmpl(10)=ngroups ! Number of groups
505  idrstmpl(11)=ngwidthref ! reference for group widths
506  idrstmpl(12)=nbitsgwidth ! num bits used for group widths
507  idrstmpl(13)=nglenref ! Reference for group lengths
508  idrstmpl(14)=1 ! length increment for group lengths
509  idrstmpl(15)=nglenlast ! True length of last group
510  idrstmpl(16)=nbitsglen ! num bits used for group lengths
511  if (idrsnum.eq.3) then
512  idrstmpl(18)=nbitsd/8 ! num bits used for extra spatial
513  ! differencing values
514  endif
515 
516  return
517  end
g2_sbytesc
subroutine g2_sbytesc(OUT, IN, ISKIP, NBYTE, NSKIP, N)
This subrountine is to put arbitrary size values into a packed bit string, taking the low order bits ...
Definition: g2_gbytesc.f:115
rdieee
subroutine rdieee(rieee, a, num)
This subroutine reads a list of real values in 32-bit IEEE floating point format.
Definition: rdieee.f:16
intmath::i1log2
Definition: intmath.f:33
g2_sbytec
subroutine g2_sbytec(OUT, IN, ISKIP, NBYTE)
This subrountine is to put arbitrary size values into a packed bit string, taking the low order bits ...
Definition: g2_gbytesc.f:39
mkieee
subroutine mkieee(a, rieee, num)
This subroutine stores a list of real values in 32-bit IEEE floating point format.
Definition: mkieee.f:17
misspack
subroutine misspack(fld, ndpts, idrsnum, idrstmpl, cpack, lcpack)
This subroutine packs up a data field using a complex packing algorithm as defined in the GRIB2 docum...
Definition: misspack.f:44
intmath
This module defines integer math functions used by other programs.
Definition: intmath.f:22
pack_gp
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:155