-
Notifications
You must be signed in to change notification settings - Fork 17
/
Copy pathcompack.F90
1301 lines (1258 loc) · 40.5 KB
/
compack.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!> @file
!> @brief Pack/unpack a data field with complex packing with or without
!> spatial differences defined in [Data Representation
!> Template5.2](https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_temp5-2.shtml)
!> and [Data Representation Template
!> 5.3](https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_temp5-3.shtml).
!> @author Stephen Gilbert @date 2000-06-21
!> @author Ed Hartnett @date Mar 5, 2024
!> Pack up a data field using a complex packing algorithm.
!>
!> This subroutine supports GRIB2 complex packing templates with or
!> without spatial differences, Data Representation Templates (DRT)
!> [GRIB2 - DATA REPRESENTATION TEMPLATE 5.2 - Grid point data -
!> complex
!> packing](https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_temp5-2.shtml)
!> and [GRIB2 - DATA REPRESENTATION TEMPLATE 5.3 - Grid point data -
!> complex packing and spatial
!> differencing](https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_temp5-3.shtml).
!>
!> It also fills in GRIB2 Data Representation Template 5.2 or 5.3 with
!> the appropriate values.
!>
!> @param[in] fld The data values to pack.
!> @param[in] ndpts The number of data values in array fld.
!> @param[in] idrsnum Data Representation Template number. Must equal
!> 2 or 3.
!> @param[inout] idrstmpl Contains the array of values for Data
!> Representation Template 5.2 or 5.3
!> - 1 Reference value - ignored on input
!> - 2 Binary Scale Factor
!> - 3 Decimal Scale Factor
!> - 7 Missing value management
!> - 8 Primary missing value
!> - 9 Secondary missing value
!> - 17 Order of Spatial Differencing (1 or 2)
!> @param[out] cpack The packed data field (character*1 array).
!> @param[out] lcpack Length of packed field cpack. -1 is returned if
!> idrstmpl(7) is not set correctly.
!>
!> @author Stephen Gilbert @date 2004-08-27
subroutine cmplxpack(fld,ndpts,idrsnum,idrstmpl,cpack,lcpack)
integer,intent(in) :: ndpts,idrsnum
real,intent(in) :: fld(ndpts)
character(len=1),intent(out) :: cpack(*)
integer,intent(inout) :: idrstmpl(*)
integer,intent(out) :: lcpack
if ( idrstmpl(7) .eq. 0 ) then ! No internal missing values
call compack(fld,ndpts,idrsnum,idrstmpl,cpack,lcpack)
elseif ( idrstmpl(7).eq.1 .OR. idrstmpl(7).eq.2) then
call misspack(fld,ndpts,idrsnum,idrstmpl,cpack,lcpack)
else
print *,'cmplxpack: Do not recognize Missing value option.'
lcpack=-1
endif
return
end subroutine cmplxpack
!> Pack a data field with complex packing with or without spatial
!> differences, defined in [Data Representation Template
!> 5.2](https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_temp5-2.shtml)
!> and [Data Representation Template
!> 5.3](https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_temp5-3.shtml).
!>
!> It also fills in Data Representation Template 5.2 or 5.3 with
!> the appropriate values.
!>
!> @param[in] fld The data values to pack.
!> @param[in] ndpts The number of data values in array fld.
!> @param[in] idrsnum Data Representation Template number - must equal
!> 2 or 3.
!> @param[inout] idrstmpl The array of values for Data
!> Representation Template 5.2 or 5.3.
!> - 1 Reference value - ignored on input
!> - 2 Binary Scale Factor
!> - 3 Decimal Scale Factor
!> - 7 Missing value management
!> - 8 Primary missing value
!> - 9 Secondary missing value
!> - 17 Order of Spatial Differencing (1 or 2)
!> @param[out] cpack The packed data field (character*1 array).
!> @param[out] lcpack length of packed field cpack.
!>
!> @author Stephen Gilbert @date 2000-06-21
subroutine compack(fld,ndpts,idrsnum,idrstmpl,cpack,lcpack)
use intmath
implicit none
integer,intent(in) :: ndpts,idrsnum
real,intent(in) :: fld(ndpts)
character(len=1),intent(out) :: cpack(*)
integer,intent(inout) :: idrstmpl(*)
integer,intent(out) :: lcpack
real :: bscale,dscale
integer :: j,iofst,imin,ival1,ival2,minsd,nbitsd,n
integer :: igmax,nbitsgref,left,iwmax,i,ilmax,kk,ij
integer :: ngwidthref,nbitsgwidth,nglenref,nglenlast
integer :: maxorig,nbitorig,isd,ngroups,itemp,minpk
integer :: inc,maxgrps,missopt,miss1,miss2,lg
integer :: ibit,jbit,kbit,novref,lbitref,ier,ng,imax
integer :: nbitsglen
real(4) :: ref,rmin4
real(8) :: rmin,rmax
integer(4) :: iref
integer,allocatable :: ifld(:)
integer,allocatable :: jmin(:),jmax(:),lbit(:)
integer,parameter :: zero=0
integer,allocatable :: gref(:),gwidth(:),glen(:)
integer :: glength,grpwidth
logical :: simple_alg
simple_alg = .false.
bscale=2.0**real(-idrstmpl(2))
dscale=10.0**real(idrstmpl(3))
!
! Find max and min values in the data
!
if(ndpts>0) then
rmax=fld(1)
rmin=fld(1)
else
rmax=1.0
rmin=1.0
endif
do j=2,ndpts
if (fld(j).gt.rmax) rmax=fld(j)
if (fld(j).lt.rmin) rmin=fld(j)
enddo
!
! If max and min values are not equal, pack up field.
! If they are equal, we have a constant field, and the reference
! value (rmin) is the value for each point in the field and
! set nbits to 0.
!
multival: if (rmin.ne.rmax) then
iofst=0
allocate(ifld(ndpts))
allocate(gref(ndpts))
allocate(gwidth(ndpts))
allocate(glen(ndpts))
!
! Scale original data
!
if (idrstmpl(2).eq.0) then ! No binary scaling
imin=nint(rmin*dscale)
!imax=nint(rmax*dscale)
rmin=real(imin)
do j=1,ndpts
ifld(j)=max(0,nint(fld(j)*dscale)-imin)
enddo
else ! Use binary scaling factor
rmin=rmin*dscale
!rmax=rmax*dscale
do j=1,ndpts
ifld(j)=max(0,nint(((fld(j)*dscale)-rmin)*bscale))
enddo
endif
!
! Calculate Spatial differences, if using DRS Template 5.3
!
alg3: if (idrsnum.eq.3) then ! spatial differences
if (idrstmpl(17).ne.1.and.idrstmpl(17).ne.2) idrstmpl(17)=2
if (idrstmpl(17).eq.1) then ! first order
ival1=ifld(1)
if(ival1<0) then
print *,'G2: negative ival1',ival1
stop 101
endif
do j=ndpts,2,-1
ifld(j)=ifld(j)-ifld(j-1)
enddo
ifld(1)=0
elseif (idrstmpl(17).eq.2) then ! second order
ival1=ifld(1)
ival2=ifld(2)
if(ival1<0) then
print *,'G2: negative ival1',ival1
stop 11
endif
if(ival2<0) then
print *,'G2: negative ival2',ival2
stop 12
endif
do j=ndpts,3,-1
ifld(j)=ifld(j)-(2*ifld(j-1))+ifld(j-2)
enddo
ifld(1)=0
ifld(2)=0
endif
!
! subtract min value from spatial diff field
!
isd=idrstmpl(17)+1
minsd=minval(ifld(isd:ndpts))
do j=isd,ndpts
ifld(j)=ifld(j)-minsd
enddo
!
! find num of bits need to store minsd and add 1 extra bit
! to indicate sign
!
nbitsd=i1log2(abs(minsd))+1
!
! find num of bits need to store ifld(1) ( and ifld(2)
! if using 2nd order differencing )
!
maxorig=ival1
if (idrstmpl(17).eq.2.and.ival2.gt.ival1) maxorig=ival2
nbitorig=i1log2(maxorig)+1
if (nbitorig.gt.nbitsd) nbitsd=nbitorig
! increase number of bits to even multiple of 8 ( octet )
if (mod(nbitsd,8).ne.0) nbitsd=nbitsd+(8-mod(nbitsd,8))
!
! Store extra spatial differencing info into the packed
! data section.
!
if (nbitsd.ne.0) then
! pack first original value
if (ival1.ge.0) then
call g2_sbytec(cpack,ival1,iofst,nbitsd)
iofst=iofst+nbitsd
else
call g2_sbytec(cpack,1,iofst,1)
iofst=iofst+1
call g2_sbytec(cpack,iabs(ival1),iofst,nbitsd-1)
iofst=iofst+nbitsd-1
endif
if (idrstmpl(17).eq.2) then
! pack second original value
if (ival2.ge.0) then
call g2_sbytec(cpack,ival2,iofst,nbitsd)
iofst=iofst+nbitsd
else
call g2_sbytec(cpack,1,iofst,1)
iofst=iofst+1
call g2_sbytec(cpack,iabs(ival2),iofst,nbitsd-1)
iofst=iofst+nbitsd-1
endif
endif
! pack overall min of spatial differences
if (minsd.ge.0) then
call g2_sbytec(cpack,minsd,iofst,nbitsd)
iofst=iofst+nbitsd
else
call g2_sbytec(cpack,1,iofst,1)
iofst=iofst+1
call g2_sbytec(cpack,iabs(minsd),iofst,nbitsd-1)
iofst=iofst+nbitsd-1
endif
endif
!print *,'SDp ',ival1,ival2,minsd,nbitsd
endif alg3 ! end of spatial diff section
!
! Determine Groups to be used.
!
simplealg: if ( simple_alg ) then
! set group length to 10 : calculate number of groups
! and length of last group
print *,'G2: use simple algorithm'
ngroups=ndpts/10
glen(1:ngroups)=10
itemp=mod(ndpts,10)
if (itemp.ne.0) then
ngroups=ngroups+1
glen(ngroups)=itemp
endif
else
! Use Dr. Glahn's algorithm for determining grouping.
!
minpk=10
inc=1
maxgrps=((ndpts+minpk-1)/minpk)
allocate(jmin(maxgrps))
allocate(jmax(maxgrps))
allocate(lbit(maxgrps))
missopt=0
call pack_gp(ifld,ndpts,missopt,minpk,inc,miss1,miss2, &
jmin,jmax,lbit,glen,maxgrps,ngroups,ibit,jbit, &
kbit,novref,lbitref,ier)
if(ier/=0) then
! Dr. Glahn's algorithm failed; use simple packing method instead.
1099 format('G2: fall back to simple algorithm (glahn ier=',I0,')')
print 1099,ier
ngroups=ndpts/10
glen(1:ngroups)=10
itemp=mod(ndpts,10)
if (itemp.ne.0) then
ngroups=ngroups+1
glen(ngroups)=itemp
endif
elseif(ngroups<1) then
! Dr. Glahn's algorithm failed; use simple packing method instead.
print *,'Glahn algorithm failed; use simple packing'
ngroups=ndpts/10
glen(1:ngroups)=10
itemp=mod(ndpts,10)
if (itemp.ne.0) then
ngroups=ngroups+1
glen(ngroups)=itemp
endif
else
!print *,'SAGier = ',ier,ibit,jbit,kbit,novref,lbitref
do ng=1,ngroups
glen(ng)=glen(ng)+novref
enddo
deallocate(jmin)
deallocate(jmax)
deallocate(lbit)
endif
endif simplealg
!
! For each group, find the group's reference value
! and the number of bits needed to hold the remaining values
!
n=1
do ng=1,ngroups
! find max and min values of group
gref(ng)=ifld(n)
imax=ifld(n)
j=n+1
do lg=2,glen(ng)
if (ifld(j).lt.gref(ng)) gref(ng)=ifld(j)
if (ifld(j).gt.imax) imax=ifld(j)
j=j+1
enddo
! calc num of bits needed to hold data
if ( gref(ng).ne.imax ) then
gwidth(ng)=i1log2(imax-gref(ng))
else
gwidth(ng)=0
endif
! Subtract min from data
j=n
do lg=1,glen(ng)
ifld(j)=ifld(j)-gref(ng)
j=j+1
enddo
! increment fld array counter
n=n+glen(ng)
enddo
!
! Find max of the group references and calc num of bits needed
! to pack each groups reference value, then
! pack up group reference values
!
!write(77,*)'GREFS: ',(gref(j),j=1,ngroups)
igmax=maxval(gref(1:ngroups))
if (igmax.ne.0) then
nbitsgref=i1log2(igmax)
call g2_sbytesc(cpack,gref,iofst,nbitsgref,0,ngroups)
itemp=nbitsgref*ngroups
iofst=iofst+itemp
if (mod(itemp,8).ne.0) then
left=8-mod(itemp,8)
call g2_sbytec(cpack,zero,iofst,left)
iofst=iofst+left
endif
else
nbitsgref=0
endif
!
! Find max/min of the group widths and calc num of bits needed
! to pack each groups width value, then
! pack up group width values
!
!write(77,*)'GWIDTHS: ',(gwidth(j),j=1,ngroups)
iwmax=maxval(gwidth(1:ngroups))
ngwidthref=minval(gwidth(1:ngroups))
if (iwmax.ne.ngwidthref) then
nbitsgwidth=i1log2(iwmax-ngwidthref)
do i=1,ngroups
gwidth(i)=gwidth(i)-ngwidthref
if(gwidth(i)<0) then
write(0,*) 'i,gw,ngw=',i,gwidth(i),ngwidthref
stop 9
endif
enddo
call g2_sbytesc(cpack,gwidth,iofst,nbitsgwidth,0,ngroups)
itemp=nbitsgwidth*ngroups
iofst=iofst+itemp
! Pad last octet with Zeros, if necessary,
if (mod(itemp,8).ne.0) then
left=8-mod(itemp,8)
call g2_sbytec(cpack,zero,iofst,left)
iofst=iofst+left
endif
else
nbitsgwidth=0
gwidth(1:ngroups)=0
endif
!
! Find max/min of the group lengths and calc num of bits needed
! to pack each groups length value, then
! pack up group length values
!
!write(77,*)'GLENS: ',(glen(j),j=1,ngroups)
ilmax=maxval(glen(1:ngroups-1))
nglenref=minval(glen(1:ngroups-1))
nglenlast=glen(ngroups)
if (ilmax.ne.nglenref) then
nbitsglen=i1log2(ilmax-nglenref)
do i=1,ngroups-1
glen(i)=glen(i)-nglenref
if(glen(i)<0) then
write(0,*) 'i,glen(i) = ',i,glen(i)
stop 23
endif
enddo
call g2_sbytesc(cpack,glen,iofst,nbitsglen,0,ngroups)
itemp=nbitsglen*ngroups
iofst=iofst+itemp
! Pad last octet with Zeros, if necessary,
if (mod(itemp,8).ne.0) then
left=8-mod(itemp,8)
call g2_sbytec(cpack,zero,iofst,left)
iofst=iofst+left
endif
else
nbitsglen=0
glen(1:ngroups)=0
endif
!
! For each group, pack data values
!
!write(77,*)'IFLDS: ',(ifld(j),j=1,ndpts)
n=1
ij=0
do ng=1,ngroups
glength=glen(ng)+nglenref
if (ng.eq.ngroups ) glength=nglenlast
grpwidth=gwidth(ng)+ngwidthref
!write(77,*)'NGP ',ng,grpwidth,glength,gref(ng)
if ( grpwidth.ne.0 ) then
call g2_sbytesc(cpack,ifld(n),iofst,grpwidth,0,glength)
iofst=iofst+(grpwidth*glength)
endif
do kk=1,glength
ij=ij+1
!write(77,*)'SAG ',ij,fld(ij),ifld(ij),gref(ng),bscale,rmin,dscale
enddo
n=n+glength
enddo
! Pad last octet with Zeros, if necessary,
if (mod(iofst,8).ne.0) then
left=8-mod(iofst,8)
call g2_sbytec(cpack,zero,iofst,left)
iofst=iofst+left
endif
lcpack=iofst/8
!
if ( allocated(ifld) ) deallocate(ifld)
if ( allocated(gref) ) deallocate(gref)
if ( allocated(gwidth) ) deallocate(gwidth)
if ( allocated(glen) ) deallocate(glen)
else ! Constant field ( max = min )
lcpack=0
nbitsgref=0
ngroups=0
ngwidthref=0
nbitsgwidth=0
nglenref=0
nglenlast=0
nbitsglen=0
nbitsd=0
endif multival
!
! Fill in ref value and number of bits in Template 5.2
!
rmin4 = real(rmin, 4)
call mkieee(rmin4,ref,1) ! ensure reference value is IEEE format
! call g2_gbytec(ref,idrstmpl(1),0,32)
iref=transfer(ref,iref)
idrstmpl(1)=iref
idrstmpl(4)=nbitsgref
idrstmpl(5)=0 ! original data were reals
idrstmpl(6)=1 ! general group splitting
idrstmpl(7)=0 ! No internal missing values
idrstmpl(8)=0 ! Primary missing value
idrstmpl(9)=0 ! secondary missing value
idrstmpl(10)=ngroups ! Number of groups
idrstmpl(11)=ngwidthref ! reference for group widths
idrstmpl(12)=nbitsgwidth ! num bits used for group widths
idrstmpl(13)=nglenref ! Reference for group lengths
idrstmpl(14)=1 ! length increment for group lengths
idrstmpl(15)=nglenlast ! True length of last group
idrstmpl(16)=nbitsglen ! num bits used for group lengths
if (idrsnum.eq.3) then
idrstmpl(18)=nbitsd/8 ! num bits used for extra spatial
! differencing values
endif
end subroutine compack
!> Unpack a data field that was packed using a complex packing
!> algorithm as defined in the GRIB2 documention.
!>
!> This subroutine Supports GRIB2 complex packing templates with or
!> without spatial differences: It supports GRIB2 complex packing
!> templates with or without spatial differences: Data Representation
!> Tables [5.2]
!> (https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_table5-2.shtml)
!> and
!> [5.3](https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_table5-3.shtml)).
!>
!> @param[in] cpack The packed data field (character*1 array).
!> @param[in] len Length of packed field cpack.
!> @param[in] lensec Length of section 7 (used for error checking).
!> @param[in] idrsnum Data Representation Template number 5.N. Must
!> equal 2 or 3.
!> @param[in] idrstmpl The array of values for Data Representation
!> Template 5.2 or 5.3.
!> @param[in] ndpts The number of data values to unpack.
!> @param[out] fld Contains the unpacked data values.
!> @param[out] ier Error return:
!> - 0 = No error.
!> - 1 = Problem - inconsistent group lengths of widths.
!>
!> @author Stephen Gilbert @date 2000-06-21
subroutine comunpack(cpack,len,lensec,idrsnum,idrstmpl,ndpts, &
fld,ier)
character(len=1),intent(in) :: cpack(len)
integer,intent(in) :: ndpts,len
integer,intent(in) :: idrstmpl(*)
real,intent(out) :: fld(ndpts)
integer,allocatable :: ifld(:),ifldmiss(:)
integer(4) :: ieee
integer,allocatable :: gref(:),gwidth(:),glen(:)
real :: ref,bscale,dscale,rmiss1,rmiss2
! real :: fldo(6045)
integer :: totBit, totLen
ier=0
!print *,'IDRSTMPL: ',(idrstmpl(j),j=1,16)
ieee = idrstmpl(1)
call rdieee(ieee,ref,1)
bscale = 2.0**real(idrstmpl(2))
dscale = 10.0**real(-idrstmpl(3))
nbitsgref = idrstmpl(4)
itype = idrstmpl(5)
ngroups = idrstmpl(10)
nbitsgwidth = idrstmpl(12)
nbitsglen = idrstmpl(16)
if (idrsnum.eq.3) then
nbitsd=idrstmpl(18)*8
endif
! Constant field
if (ngroups.eq.0) then
do j=1,ndpts
fld(j)=ref
enddo
return
endif
iofst=0
allocate(ifld(ndpts),stat=is)
!print *,'ALLOC ifld: ',is,ndpts
allocate(gref(ngroups),stat=is)
!print *,'ALLOC gref: ',is
allocate(gwidth(ngroups),stat=is)
!print *,'ALLOC gwidth: ',is
!
! Get missing values, if supplied
!
if ( idrstmpl(7).eq.1 ) then
if (itype.eq.0) then
call rdieee(idrstmpl(8),rmiss1,1)
else
rmiss1=real(idrstmpl(8))
endif
elseif ( idrstmpl(7).eq.2 ) then
if (itype.eq.0) then
call rdieee(idrstmpl(8),rmiss1,1)
call rdieee(idrstmpl(9),rmiss2,1)
else
rmiss1=real(idrstmpl(8))
rmiss2=real(idrstmpl(9))
endif
endif
!print *,'RMISSs: ',rmiss1,rmiss2,ref
!
! Extract Spatial differencing values, if using DRS Template 5.3
!
if (idrsnum.eq.3) then
if (nbitsd.ne.0) then
call g2_gbytec(cpack,ival1,iofst,nbitsd)
iofst=iofst+nbitsd
if (idrstmpl(17).eq.2) then
call g2_gbytec(cpack,ival2,iofst,nbitsd)
iofst=iofst+nbitsd
endif
call g2_gbytec(cpack,isign,iofst,1)
iofst=iofst+1
call g2_gbytec(cpack,minsd,iofst,nbitsd-1)
iofst=iofst+nbitsd-1
if (isign.eq.1) minsd=-minsd
else
ival1=0
ival2=0
minsd=0
endif
!print *,'SDu ',ival1,ival2,minsd,nbitsd
endif
!
! Extract Each Group's reference value
!
!print *,'SAG1: ',nbitsgref,ngroups,iofst
if (nbitsgref.ne.0) then
call g2_gbytesc(cpack,gref,iofst,nbitsgref,0,ngroups)
itemp=nbitsgref*ngroups
iofst=iofst+(itemp)
if (mod(itemp,8).ne.0) iofst=iofst+(8-mod(itemp,8))
else
gref(1:ngroups)=0
endif
!write(78,*)'GREFs: ',(gref(j),j=1,ngroups)
!
! Extract Each Group's bit width
!
!print *,'SAG2: ',nbitsgwidth,ngroups,iofst,idrstmpl(11)
if (nbitsgwidth.ne.0) then
call g2_gbytesc(cpack,gwidth,iofst,nbitsgwidth,0,ngroups)
itemp=nbitsgwidth*ngroups
iofst=iofst+(itemp)
if (mod(itemp,8).ne.0) iofst=iofst+(8-mod(itemp,8))
else
gwidth(1:ngroups)=0
endif
do j=1,ngroups
gwidth(j)=gwidth(j)+idrstmpl(11)
enddo
!write(78,*)'GWIDTHs: ',(gwidth(j),j=1,ngroups)
!
! Extract Each Group's length (number of values in each group)
!
allocate(glen(ngroups),stat=is)
!print *,'ALLOC glen: ',is
!print *,'SAG3: ',nbitsglen,ngroups,iofst,idrstmpl(14),idrstmpl(13)
if (nbitsglen.ne.0) then
call g2_gbytesc(cpack,glen,iofst,nbitsglen,0,ngroups)
itemp=nbitsglen*ngroups
iofst=iofst+(itemp)
if (mod(itemp,8).ne.0) iofst=iofst+(8-mod(itemp,8))
else
glen(1:ngroups)=0
endif
do j=1,ngroups
glen(j)=(glen(j)*idrstmpl(14))+idrstmpl(13)
enddo
glen(ngroups)=idrstmpl(15)
!write(78,*)'GLENs: ',(glen(j),j=1,ngroups)
!print *,'GLENsum: ',sum(glen)
!
! Test to see if the group widths and lengths are consistent with number of
! values, and length of section 7.
!
totBit = 0
totLen = 0
do j=1,ngroups
totBit = totBit + (gwidth(j)*glen(j));
totLen = totLen + glen(j);
enddo
if (totLen .NE. ndpts) then
ier=1
return
endif
if ( (totBit/8) .GT. lensec) then
ier=1
return
endif
!
! For each group, unpack data values
!
if ( idrstmpl(7).eq.0 ) then ! no missing values
n=1
do j=1,ngroups
!write(78,*)'NGP ',j,gwidth(j),glen(j),gref(j)
if (gwidth(j).ne.0) then
call g2_gbytesc(cpack,ifld(n),iofst,gwidth(j),0,glen(j))
do k=1,glen(j)
ifld(n)=ifld(n)+gref(j)
n=n+1
enddo
else
ifld(n:n+glen(j)-1)=gref(j)
n=n+glen(j)
endif
iofst=iofst+(gwidth(j)*glen(j))
enddo
elseif ( idrstmpl(7).eq.1.OR.idrstmpl(7).eq.2 ) then
! missing values included
allocate(ifldmiss(ndpts))
!ifldmiss=0
n=1
non=1
do j=1,ngroups
!print *,'SAGNGP ',j,gwidth(j),glen(j),gref(j)
if (gwidth(j).ne.0) then
msng1=(2**gwidth(j))-1
msng2=msng1-1
call g2_gbytesc(cpack,ifld(n),iofst,gwidth(j),0,glen(j))
iofst=iofst+(gwidth(j)*glen(j))
do k=1,glen(j)
if (ifld(n).eq.msng1) then
ifldmiss(n)=1
elseif (idrstmpl(7).eq.2.AND.ifld(n).eq.msng2) then
ifldmiss(n)=2
else
ifldmiss(n)=0
ifld(non)=ifld(n)+gref(j)
non=non+1
endif
n=n+1
enddo
else
msng1=(2**nbitsgref)-1
msng2=msng1-1
if (gref(j).eq.msng1) then
ifldmiss(n:n+glen(j)-1)=1
!ifld(n:n+glen(j)-1)=0
elseif (idrstmpl(7).eq.2.AND.gref(j).eq.msng2) then
ifldmiss(n:n+glen(j)-1)=2
!ifld(n:n+glen(j)-1)=0
else
ifldmiss(n:n+glen(j)-1)=0
ifld(non:non+glen(j)-1)=gref(j)
non=non+glen(j)
endif
n=n+glen(j)
endif
enddo
endif
!write(78,*)'IFLDs: ',(ifld(j),j=1,ndpts)
if ( allocated(gref) ) deallocate(gref)
if ( allocated(gwidth) ) deallocate(gwidth)
if ( allocated(glen) ) deallocate(glen)
!
! If using spatial differences, add overall min value, and
! sum up recursively
!
if (idrsnum.eq.3) then ! spatial differencing
if (idrstmpl(17).eq.1) then ! first order
ifld(1)=ival1
if ( idrstmpl(7).eq.0 ) then ! no missing values
itemp=ndpts
else
itemp=non-1
endif
do n=2,itemp
ifld(n)=ifld(n)+minsd
ifld(n)=ifld(n)+ifld(n-1)
enddo
elseif (idrstmpl(17).eq.2) then ! second order
ifld(1)=ival1
ifld(2)=ival2
if ( idrstmpl(7).eq.0 ) then ! no missing values
itemp=ndpts
else
itemp=non-1
endif
do n=3,itemp
ifld(n)=ifld(n)+minsd
ifld(n)=ifld(n)+(2*ifld(n-1))-ifld(n-2)
enddo
endif
!write(78,*)'IFLDs: ',(ifld(j),j=1,ndpts)
endif
!
! Scale data back to original form
!
!print *,'SAGT: ',ref,bscale,dscale
if ( idrstmpl(7).eq.0 ) then ! no missing values
do n=1,ndpts
fld(n)=((real(ifld(n))*bscale)+ref)*dscale
!write(78,*)'SAG ',n,fld(n),ifld(n),bscale,ref,1./dscale
enddo
elseif ( idrstmpl(7).eq.1.OR.idrstmpl(7).eq.2 ) then
! missing values included
non=1
do n=1,ndpts
if ( ifldmiss(n).eq.0 ) then
fld(n)=((real(ifld(non))*bscale)+ref)*dscale
!print *,'SAG ',n,fld(n),ifld(non),bscale,ref,dscale
non=non+1
elseif ( ifldmiss(n).eq.1 ) then
fld(n)=rmiss1
elseif ( ifldmiss(n).eq.2 ) then
fld(n)=rmiss2
endif
enddo
if ( allocated(ifldmiss) ) deallocate(ifldmiss)
endif
if ( allocated(ifld) ) deallocate(ifld)
!open(10,form='unformatted',recl=24180,access='direct')
!read(10,rec=1) (fldo(k),k=1,6045)
!do i =1,6045
! print *,i,fldo(i),fld(i),fldo(i)-fld(i)
!enddo
return
end subroutine comunpack
!> Pack up a data field using a GRIB2 algorithm with missing value
!> management.
!>
!> This subroutine packs up a data field using a complex packing
!> algorithm as defined in the GRIB2 documention. It supports GRIB2
!> complex packing templates with or without spatial differences
!> (i.e. Data Representation Tables [5.2]
!> (https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_table5-2.shtml)
!> and
!> [5.3](https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_table5-3.shtml)).
!>
!> This subroutine also fills in GRIB2 Data Representation Template
!> 5.2 or 5.3 with the appropriate values. It assumes that Missing
!> Value Management is being used and that 1 or 2 missing values
!> appear in the data.
!>
!> @param[in] fld Contains the data values to pack.
!> @param[in] ndpts The number of data values in array fld.
!> @param[in] idrsnum Data Representation Template number. Must
!> equal 2 or 3.
!> @param[inout] idrstmpl Contains the array of values for Data
!> Representation Template 5.2 or 5.3
!> - idrstmpl(1) Reference value - ignored on input set by compack routine.
!> - idrstmpl(2) Binary Scale Factor.
!> - idrstmpl(3) Decimal Scale Factor.
!> - idrstmpl(4) number of bits for each data value - ignored on input.
!> - idrstmpl(5) Original field type, currently ignored on input, set = 0 on
!> !output, Data values assumed to be reals.
!> - idrstmpl(6) = 0 use lossless compression or = 1 use lossy compression.
!> - idrstmpl(7) Missing value management.
!> - idrstmpl(8) Primary missing value.
!> - idrstmpl(9) Secondary missing value.
!> - idrstmpl(17) Order of Spatial Differencing (1 or 2).
!> @param[out] cpack The packed data field (character*1 array).
!> @param[out] lcpack length of packed field cpack. -1 is returned if
!> idrstmpl(7) is not set correctly.
!>
!> @author Stephen Gilbert @date 2000-06-21
subroutine misspack(fld,ndpts,idrsnum,idrstmpl,cpack,lcpack)
use intmath
integer,intent(in) :: ndpts,idrsnum
real,intent(in) :: fld(ndpts)
character(len=1),intent(out) :: cpack(*)
integer,intent(inout) :: idrstmpl(*)
integer,intent(out) :: lcpack
real(4) :: ref, rmin4
integer(4) :: iref
integer,allocatable :: ifld(:),ifldmiss(:),jfld(:)
integer,allocatable :: jmin(:),jmax(:),lbit(:)
integer,parameter :: zero=0
integer,allocatable :: gref(:),gwidth(:),glen(:)
integer :: glength,grpwidth
logical :: simple_alg
logical :: have_rmin
simple_alg = .false.
have_rmin = .false.
bscale=2.0**real(-idrstmpl(2))
dscale=10.0**real(idrstmpl(3))
missopt=idrstmpl(7)
if ( missopt.ne.1 .AND. missopt.ne.2 ) then
print *,'misspack: Unrecognized option.'
lcpack=-1
return
else ! Get missing values
call rdieee(idrstmpl(8),rmissp,1)
if (missopt.eq.2) call rdieee(idrstmpl(9),rmisss,1)
endif
! Find min value of non-missing values in the data,
! AND set up missing value mapping of the field.
allocate(ifldmiss(ndpts))
if ( missopt .eq. 1 ) then ! Primary missing value only
do j=1,ndpts
if (fld(j).eq.rmissp) then
ifldmiss(j)=1
else
ifldmiss(j)=0
if(have_rmin) then
if (fld(j).lt.rmin) rmin=fld(j)
else
rmin=fld(j)
have_rmin=.true.
endif
endif
enddo
if(.not.have_rmin) rmin=rmissp
endif
if ( missopt .eq. 2 ) then ! Primary and secondary missing values
do j=1,ndpts
if (fld(j).eq.rmissp) then
ifldmiss(j)=1
elseif (fld(j).eq.rmisss) then
ifldmiss(j)=2
else
ifldmiss(j)=0
if(have_rmin) then
if (fld(j).lt.rmin) rmin=fld(j)
else
rmin=fld(j)
have_rmin=.true.
endif
endif
if(.not.have_rmin) rmin=rmissp
enddo
endif
! Allocate work arrays:
! Note: -ifldmiss(j),j=1,ndpts is a map of original field indicating
! which of the original data values
! are primary missing (1), sencondary missing (2) or non-missing (0).
! -jfld(j),j=1,nonmiss is a subarray of just the non-missing values from
! the original field.
iofst=0
allocate(ifld(ndpts))
allocate(jfld(ndpts))
allocate(gref(ndpts))
allocate(gwidth(ndpts))
allocate(glen(ndpts))
!
! Scale original data
!
nonmiss=0
if (idrstmpl(2).eq.0) then ! No binary scaling
imin=nint(rmin*dscale)
!imax=nint(rmax*dscale)
rmin=real(imin)
do j=1,ndpts
if (ifldmiss(j).eq.0) then
nonmiss=nonmiss+1
jfld(nonmiss)=max(0,nint(fld(j)*dscale)-imin)
endif
enddo
else ! Use binary scaling factor
rmin=rmin*dscale
do j=1,ndpts
if (ifldmiss(j).eq.0) then
nonmiss=nonmiss+1
jfld(nonmiss)=max(0,nint(((fld(j)*dscale)-rmin)*bscale))
endif
enddo
endif
!
! Calculate Spatial differences, if using DRS Template 5.3
!
if (idrsnum.eq.3) then ! spatial differences
if (idrstmpl(17).ne.1.and.idrstmpl(17).ne.2) idrstmpl(17)=2
if (idrstmpl(17).eq.1) then ! first order
if(nonmiss<1) then
ival1=1.0
else
ival1=jfld(1)
endif
do j=nonmiss,2,-1
jfld(j)=jfld(j)-jfld(j-1)
enddo
if(nonmiss>0) jfld(1)=0
elseif (idrstmpl(17).eq.2) then ! second order
if(nonmiss==1) then
ival1=jfld(1)
ival2=jfld(1)
elseif(nonmiss<1) then
ival1=1.0
ival2=1.0
else
ival1=jfld(1)
ival2=jfld(2)
endif
do j=nonmiss,3,-1
jfld(j)=jfld(j)-(2*jfld(j-1))+jfld(j-2)
enddo
if(nonmiss>=1) jfld(1)=0
if(nonmiss>=2) jfld(2)=0
endif
! subtract min value from spatial diff field
isd=idrstmpl(17)+1
minsd=minval(jfld(isd:nonmiss))
do j=isd,nonmiss
jfld(j)=jfld(j)-minsd
enddo