NCEPLIBS-g2  3.4.7
pngpack.F90
Go to the documentation of this file.
1 
4 
28 subroutine pngpack(fld, width, height, idrstmpl, cpack, lcpack)
29  implicit none
30 
31  integer, intent(in) :: width, height
32  real, intent(in) :: fld(width * height)
33  character(len = 1), intent(out) :: cpack(*)
34  integer, intent(inout) :: idrstmpl(*)
35  integer, intent(out) :: lcpack
36 
37  real(4) :: ref, rmin4
38  real(8) :: rmin, rmax
39  integer(4) :: iref
40  integer :: ifld(width * height), nbits
41  integer, parameter :: zero = 0
42  character(len = 1), allocatable :: ctemp(:)
43  real :: bscale, dscale, temp
44  integer :: imax, imin, j, maxdif, nbytes, ndpts
45 
46  interface
47  function enc_png(data, width, height, nbits, pngbuf) bind(c, name="enc_png")
48  use iso_c_binding
49  character(kind = c_char), intent(in) :: data(*)
50  integer(c_int), intent(in) :: width, height
51  integer(c_int), intent(inout) :: nbits
52  character(kind = c_char), intent(out) :: pngbuf(*)
53  integer(c_int) :: enc_png
54  end function enc_png
55  end interface
56 
57  ndpts = width * height
58  bscale = 2.0**real(-idrstmpl(2))
59  dscale = 10.0**real(idrstmpl(3))
60 
61  ! Find max and min values in the data
62  if(ndpts > 0) then
63  rmax = fld(1)
64  rmin = fld(1)
65  else
66  rmax = 1.0
67  rmin = 1.0
68  endif
69  do j = 2, ndpts
70  if (fld(j) .gt. rmax) rmax = fld(j)
71  if (fld(j) .lt. rmin) rmin = fld(j)
72  enddo
73  maxdif = nint((rmax - rmin) * dscale * bscale)
74 
75  ! If max and min values are not equal, pack up field. If they are
76  ! equal, we have a constant field, and the reference value (rmin) is
77  ! the value for each point in the field and set nbits to 0.
78  if (rmin .ne. rmax .AND. maxdif .ne. 0) then
79 
80  ! Determine which algorithm to use based on user-supplied binary
81  ! scale factor and number of bits.
82  if (idrstmpl(2) .eq. 0) then
83  ! No binary scaling and calculate minimum number of bits in
84  ! which the data will fit.
85  imin = nint(rmin * dscale)
86  imax = nint(rmax * dscale)
87  maxdif = imax - imin
88  temp = alog(real(maxdif + 1)) / alog(2.0)
89  nbits = ceiling(temp)
90  rmin = real(imin)
91  ! scale data
92  do j = 1, ndpts
93  ifld(j) = nint(fld(j) * dscale) - imin
94  enddo
95  else
96  ! Use binary scaling factor and calculate minimum number of
97  ! bits in which the data will fit.
98  rmin = rmin * dscale
99  rmax = rmax * dscale
100  maxdif = nint((rmax - rmin) * bscale)
101  temp = alog(real(maxdif + 1)) / alog(2.0)
102  nbits = ceiling(temp)
103  ! scale data
104  do j = 1, ndpts
105  ifld(j) = max(0, nint(((fld(j) * dscale) - rmin) * bscale))
106  enddo
107  endif
108 
109  ! Pack data into full octets, then do PNG encode. and calculate
110  ! the length of the packed data in bytes.
111  if (nbits .le. 8) then
112  nbits = 8
113  elseif (nbits .le. 16) then
114  nbits = 16
115  elseif (nbits .le. 24) then
116  nbits = 24
117  else
118  nbits = 32
119  endif
120  nbytes = (nbits / 8) * ndpts
121  allocate(ctemp(nbytes))
122  call g2_sbytesc(ctemp, ifld, 0, nbits, 0, ndpts)
123 
124  ! Encode data into PNG Format.
125  lcpack = enc_png(ctemp, width, height, nbits, cpack)
126  if (lcpack .le. 0) then
127  print *, 'pngpack: ERROR Encoding PNG = ', lcpack
128  endif
129  deallocate(ctemp)
130 
131  else
132  nbits = 0
133  lcpack = 0
134  endif
135 
136  ! Fill in ref value and number of bits in Template 5.0.
137  rmin4 = real(rmin, 4)
138  call mkieee(rmin4, ref, 1) ! ensure reference value is IEEE format
139  iref = transfer(ref, iref)
140  idrstmpl(1) = iref
141  idrstmpl(4) = nbits
142  idrstmpl(5) = 0 ! original data were reals
143 
144 end subroutine pngpack
int enc_png(char *data, g2int *width, g2int *height, g2int *nbits, char *pngbuf)
create png_structs to write png stream into memory.
Definition: enc_png.c:79
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
subroutine pngpack(fld, width, height, idrstmpl, cpack, lcpack)
This subroutine packs up a data field into PNG image format.
Definition: pngpack.F90:29