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