NCEPLIBS-g2  3.4.7
specpack.F90
Go to the documentation of this file.
1 
4 
21 subroutine specpack(fld,ndpts,JJ,KK,MM,idrstmpl,cpack,lcpack)
22 
23  real,intent(in) :: fld(ndpts)
24  integer,intent(in) :: ndpts,JJ,KK,MM
25  integer,intent(inout) :: idrstmpl(*)
26  character(len=1),intent(out) :: cpack(*)
27  integer,intent(out) :: lcpack
28 
29  integer :: Ts,tmplsim(5)
30  real :: bscale,dscale,unpk(ndpts),tfld(ndpts)
31  real,allocatable :: pscale(:)
32 
33  bscale = 2.0**real(-idrstmpl(2))
34  dscale = 10.0**real(idrstmpl(3))
35  nbits = idrstmpl(4)
36  js=idrstmpl(6)
37  ks=idrstmpl(7)
38  ms=idrstmpl(8)
39  ts=idrstmpl(9)
40 
41  ! Calculate Laplacian scaling factors for each possible wave number.
42  allocate(pscale(jj+mm))
43  tscale=real(idrstmpl(5))*1e-6
44  do n=js,jj+mm
45  pscale(n)=real(n*(n+1))**(tscale)
46  enddo
47 
48  ! Separate spectral coeffs into two lists; one to contain unpacked
49  ! values within the sub-spectrum Js, Ks, Ms, and the other with
50  ! values outside of the sub-spectrum to be packed.
51  inc=1
52  incu=1
53  incp=1
54  do m=0,mm
55  nm=jj ! triangular or trapezoidal
56  if (kk .eq. jj+mm) nm=jj+m ! rhombodial
57  ns=js ! triangular or trapezoidal
58  if (ks .eq. js+ms) ns=js+m ! rhombodial
59  do n=m,nm
60  if (n.le.ns .AND. m.le.ms) then ! save unpacked value
61  unpk(incu)=fld(inc) ! real part
62  unpk(incu+1)=fld(inc+1) ! imaginary part
63  inc=inc+2
64  incu=incu+2
65  else
66  ! Save value to be packed and scale Laplacian scale factor.
67  tfld(incp)=fld(inc)*pscale(n) ! real part
68  tfld(incp+1)=fld(inc+1)*pscale(n) ! imaginary part
69  inc=inc+2
70  incp=incp+2
71  endif
72  enddo
73  enddo
74 
75  deallocate(pscale)
76 
77  incu=incu-1
78  if (incu .ne. ts) then
79  print *,'specpack: Incorrect number of unpacked values ', 'given:',ts
80  print *,'specpack: Resetting idrstmpl(9) to ',incu
81  ts=incu
82  endif
83 
84  ! Add unpacked values to the packed data array in 32-bit IEEE format.
85  call mkieee(unpk,cpack,ts)
86  ipos=4*ts
87 
88  ! Scale and pack the rest of the coefficients.
89  tmplsim(2)=idrstmpl(2)
90  tmplsim(3)=idrstmpl(3)
91  tmplsim(4)=idrstmpl(4)
92  call simpack(tfld,ndpts-ts,tmplsim,cpack(ipos+1),lcpack)
93  lcpack=lcpack+ipos
94 
95  ! Fill in Template 5.51.
96  idrstmpl(1)=tmplsim(1)
97  idrstmpl(2)=tmplsim(2)
98  idrstmpl(3)=tmplsim(3)
99  idrstmpl(4)=tmplsim(4)
100  idrstmpl(9)=ts
101  idrstmpl(10)=1 ! Unpacked spectral data is 32-bit IEEE
102 
103 end subroutine specpack
subroutine mkieee(a, rieee, num)
Copy an array of real to an array of 32-bit IEEE floating points.
Definition: mkieee.F90:14
subroutine simpack(fld, ndpts, idrstmpl, cpack, lcpack)
Pack up a data field using a simple packing algorithm as defined in the GRIB2 documention.
Definition: simpack.f:26
subroutine specpack(fld, ndpts, JJ, KK, MM, idrstmpl, cpack, lcpack)
This subroutine packs a spectral data field using the complex packing algorithm for spherical harmoni...
Definition: specpack.F90:22