NCEPLIBS-g2  3.4.9
g2bytes.F90
Go to the documentation of this file.
1 
5 
20 subroutine g2_gbytec(in, iout, iskip, nbits)
21  implicit none
22 
23  character*1, intent(in) :: in(*)
24  integer, intent(inout) :: iout(*)
25  integer, intent(in) :: iskip, nbits
26  call g2_gbytesc(in, iout, iskip, nbits, 0, 1)
27 end subroutine g2_gbytec
28 
42 subroutine g2_gbytec1(in, siout, iskip, nbits)
43  implicit none
44 
45  character*1, intent(in) :: in(*)
46  integer, intent(inout) :: siout
47  integer, intent(in) :: iskip, nbits
48  integer (kind = 4) :: iout(1)
49 
50  call g2_gbytesc(in, iout, iskip, nbits, 0, 1)
51  siout = iout(1)
52 end subroutine g2_gbytec1
53 
63 !iteration.
67 subroutine g2_gbytescr(in, rout, iskip, nbits, nskip, n)
68  implicit none
69  character*1, intent(in) :: in(*)
70  real (kind = 4), intent(out) :: rout(*)
71  integer, intent(in) :: iskip, nbits, nskip, n
72  integer (kind = 4) :: iout(n)
73 
74  ! Unpack into integer array.
75  call g2_gbytesc(in, iout, iskip, nbits, nskip, n)
76 
77  ! Transfer to real array.
78  rout(1:n) = transfer(iout, rout(1:n), n)
79 end subroutine g2_gbytescr
80 
93 subroutine g2_gbytesc(in, iout, iskip, nbits, nskip, n)
94  implicit none
95 
96  character*1, intent(in) :: in(*)
97  integer, intent(out) :: iout(*)
98  integer, intent(in) :: iskip, nbits, nskip, n
99  integer :: tbit, bitcnt
100  integer, parameter :: ones(8) = (/ 1, 3, 7, 15, 31, 63, 127, 255 /)
101 
102  integer :: nbit, i, index, ibit, itmp
103  integer, external :: g2_mova2i
104 
105  ! nbit is the start position of the field in bits
106  nbit = iskip
107  do i = 1, n
108  bitcnt = nbits
109  index = nbit / 8 + 1
110  ibit = mod(nbit, 8)
111  nbit = nbit + nbits + nskip
112 
113  ! first byte
114  tbit = min(bitcnt, 8 - ibit)
115  itmp = iand(g2_mova2i(in(index)), ones(8 - ibit))
116  if (tbit .ne. 8 - ibit) itmp = ishft(itmp, tbit - 8 + ibit)
117  index = index + 1
118  bitcnt = bitcnt - tbit
119 
120  ! now transfer whole bytes
121  do while (bitcnt .ge. 8)
122  itmp = ior(ishft(itmp,8), g2_mova2i(in(index)))
123  bitcnt = bitcnt - 8
124  index = index + 1
125  enddo
126 
127  ! get data from last byte
128  if (bitcnt .gt. 0) then
129  itmp = ior(ishft(itmp, bitcnt), iand(ishft(g2_mova2i(in(index)), &
130  - (8 - bitcnt)), ones(bitcnt)))
131  endif
132 
133  iout(i) = itmp
134  enddo
135 
136 end subroutine g2_gbytesc
137 
151 subroutine g2_gbytec8(in, iout, iskip, nbits)
152  implicit none
153 
154  character*1, intent(in) :: in(*)
155  integer (kind = 8), intent(inout) :: iout(*)
156  integer, intent(in) :: iskip, nbits
157  call g2_gbytesc8(in, iout, iskip, nbits, 0, 1)
158 end subroutine g2_gbytec8
159 
172 subroutine g2_gbytesc8(in, iout, iskip, nbits, nskip, n)
173  implicit none
174 
175  character*1, intent(in) :: in(*)
176  integer (kind = 8), intent(out) :: iout(*)
177  integer, intent(in) :: iskip, nbits, nskip, n
178  integer :: tbit, bitcnt
179  integer, parameter :: ones(8) = (/ 1, 3, 7, 15, 31, 63, 127, 255 /)
180 
181  integer :: nbit, i, index, ibit, itmp
182  integer (kind = 8) :: itmp8, itmp8_2, itmp8_3
183  integer, external :: g2_mova2i
184  integer (kind = 8), external :: g2_mova2i8
185 
186  ! nbit is the start position of the field in bits
187  nbit = iskip
188  do i = 1, n
189  bitcnt = nbits
190  index = nbit / 8 + 1
191  ibit = mod(nbit, 8)
192  nbit = nbit + nbits + nskip
193 
194  ! first byte
195  tbit = min(bitcnt, 8 - ibit)
196  itmp8 = iand(g2_mova2i8(in(index)), int(ones(8 - ibit), kind = 8))
197  itmp = iand(g2_mova2i(in(index)), ones(8 - ibit))
198  if (tbit .ne. 8 - ibit) itmp = ishft(itmp, tbit - 8 + ibit)
199  if (tbit .ne. 8 - ibit) itmp8 = ishft(itmp8, tbit - 8 + ibit)
200  index = index + 1
201  bitcnt = bitcnt - tbit
202 
203  ! now transfer whole bytes
204  do while (bitcnt .ge. 8)
205  itmp = ior(ishft(itmp,8), g2_mova2i(in(index)))
206  itmp8 = ior(ishft(itmp8,8), g2_mova2i8(in(index)))
207  bitcnt = bitcnt - 8
208  index = index + 1
209  enddo
210 
211  ! get data from last byte
212  if (bitcnt .gt. 0) then
213  itmp = ior(ishft(itmp, bitcnt), iand(ishft(g2_mova2i(in(index)), - (8 - bitcnt)), ones(bitcnt)))
214  itmp8_2 = ishft(g2_mova2i8(in(index)), int(-(8 - bitcnt), kind(8)))
215  itmp8_3 = int(ones(bitcnt), kind(8))
216  itmp8 = ior(ishft(itmp8, bitcnt), iand(itmp8_2, itmp8_3))
217  endif
218 
219  iout(i) = itmp8
220  enddo
221 
222 end subroutine g2_gbytesc8
223 
237 subroutine g2_sbytec(out, in, iskip, nbits)
238  implicit none
239 
240  character*1, intent(inout) :: out(*)
241  integer, intent(in) :: in(*)
242  integer, intent(in) :: iskip, nbits
243  call g2_sbytesc(out, in, iskip, nbits, 0, 1)
244 end subroutine g2_sbytec
245 
259 subroutine g2_sbytec1(out, in, iskip, nbits)
260  implicit none
261 
262  character*1, intent(inout) :: out(*)
263  integer, intent(in) :: in
264  integer, intent(in) :: iskip, nbits
265  integer :: ain(1)
266  ain(1) = in
267  call g2_sbytesc(out, ain, iskip, nbits, 0, 1)
268 end subroutine g2_sbytec1
269 
281 subroutine g2_sbytescr(out, rin, iskip, nbits, nskip, n)
282  implicit none
283  character*1, intent(out) :: out(*)
284  real, intent(in) :: rin(n)
285  integer, intent(in) :: iskip, nbits, nskip, n
286  integer :: in(n)
287 
288  ! Transfer real array to integer array.
289  in(1:n) = transfer(rin, in(1:n), n)
290 
291  ! Pack into character array.
292  call g2_sbytesc(out, in, iskip, nbits, nskip, n)
293 end subroutine g2_sbytescr
294 
307 subroutine g2_sbytesc(out, in, iskip, nbits, nskip, n)
308  implicit none
309 
310  character*1, intent(out) :: out(*)
311  integer, intent(in) :: in(n)
312  integer, intent(in) :: iskip, nbits, nskip, n
313 
314  integer :: bitcnt, tbit
315  integer, parameter :: ones(8)=(/ 1, 3, 7, 15, 31, 63, 127, 255/)
316  integer :: nbit, i, itmp, index, ibit, imask, itmp2, itmp3
317  integer, external :: g2_mova2i
318 
319  ! number bits from zero to ...
320  ! nbit is the last bit of the field to be filled
321  nbit = iskip + nbits - 1
322  !print *, 'nbit', nbit, 'nbits ', nbits, 'nskip', nskip, 'n', n
323  do i = 1, n
324  itmp = in(i)
325  bitcnt = nbits
326  index = nbit / 8 + 1
327  ibit = mod(nbit, 8)
328  nbit = nbit + nbits + nskip
329  !print *, 'i', i, 'itmp', itmp, 'bitcnt', bitcnt, 'index', index, 'ibit', ibit, 'nbit', nbit
330 
331  ! make byte aligned
332  if (ibit .ne. 7) then
333  tbit = min(bitcnt, ibit + 1)
334  imask = ishft(ones(tbit), 7 - ibit)
335  itmp2 = iand(ishft(itmp, 7 - ibit),imask)
336  itmp3 = iand(g2_mova2i(out(index)), 255 - imask)
337  out(index) = char(ior(itmp2, itmp3))
338  bitcnt = bitcnt - tbit
339  itmp = ishft(itmp, -tbit)
340  index = index - 1
341  endif
342 
343  ! now byte aligned
344 
345  ! do by bytes
346  do while (bitcnt .ge. 8)
347  out(index) = char(iand(itmp, 255))
348  !print '(z2.2, x, z2.2, x, z2.2)', out(index), itmp, iand(itmp, 255)
349  itmp = ishft(itmp, -8)
350  bitcnt = bitcnt - 8
351  index = index - 1
352  enddo
353 
354  ! Do left over bits.
355  if (bitcnt .gt. 0) then
356  itmp2 = iand(itmp, ones(bitcnt))
357  !print '(z2.2, x, z2.2)', ones(bitcnt), itmp2
358  itmp3 = iand(g2_mova2i(out(index)), 255 - ones(bitcnt))
359  out(index) = char(ior(itmp2, itmp3))
360  endif
361  enddo
362 
363 end subroutine g2_sbytesc
364 
379 subroutine g2_sbytec8(out, in, iskip, nbits)
380  implicit none
381 
382  character*1, intent(inout) :: out(*)
383  integer (kind = 8), intent(in) :: in(*)
384  integer, intent(in) :: iskip, nbits
385  call g2_sbytesc8(out, in, iskip, nbits, 0, 1)
386 end subroutine g2_sbytec8
387 
401 subroutine g2_sbytesc8(out, in, iskip, nbits, nskip, n)
402  implicit none
403 
404  character*1, intent(out) :: out(*)
405  integer (kind = 8), intent(in) :: in(n)
406  integer, intent(in) :: iskip, nbits, nskip, n
407 
408  integer :: bitcnt, tbit
409  integer, parameter :: ones(8)=(/ 1, 3, 7, 15, 31, 63, 127, 255/)
410  integer :: nbit, i, index, ibit, imask, itmp1, itmp2, itmp3
411  integer (kind = 8) :: itmp8, itmp8_2
412  integer, external :: g2_mova2i
413 
414  ! number bits from zero to ...
415  ! nbit is the last bit of the field to be filled
416  nbit = iskip + nbits - 1
417  !print *, 'nbit', nbit, 'nbits ', nbits, 'nskip', nskip, 'n', n
418  do i = 1, n
419  itmp8 = in(i)
420  bitcnt = nbits
421  index = nbit / 8 + 1
422  ibit = mod(nbit, 8)
423  nbit = nbit + nbits + nskip
424  !print *, 'i', i, 'itmp8', itmp8, 'bitcnt', bitcnt, 'index', index, 'ibit', ibit, 'nbit', nbit
425 
426  ! make byte aligned
427  if (ibit .ne. 7) then
428  tbit = min(bitcnt, ibit + 1)
429  imask = ishft(ones(tbit), 7 - ibit)
430  itmp1 = int(ishft(itmp8, int(7 - ibit, kind(8))), kind(4))
431  itmp2 = iand(itmp1, imask)
432  itmp3 = iand(g2_mova2i(out(index)), 255 - imask)
433  out(index) = char(ior(itmp2, itmp3))
434  bitcnt = bitcnt - tbit
435  itmp8 = ishft(itmp8, -tbit)
436  index = index - 1
437  endif
438 
439  ! now byte aligned
440 
441  ! Process a byte at a time.
442  do while (bitcnt .ge. 8)
443  !print *, bitcnt, iand(itmp8, 255_8)
444  out(index) = char(iand(itmp8, 255_8))
445  itmp8 = ishft(itmp8, -8)
446  bitcnt = bitcnt - 8
447  index = index - 1
448  enddo
449 
450  ! Take care of left over bits.
451  if (bitcnt .gt. 0) then
452  itmp8_2 = int(ones(bitcnt), kind(8))
453  itmp2 = int(iand(itmp8, itmp8_2), kind(4))
454  itmp3 = iand(g2_mova2i(out(index)), 255 - ones(bitcnt))
455  out(index) = char(ior(itmp2, itmp3))
456  endif
457  enddo
458 end subroutine g2_sbytesc8
459 
469 subroutine rdieeec(cieee, a, num)
470  implicit none
471 
472  character(len = 1), intent(in) :: cieee(*)
473  real, intent(out) :: a(num)
474  integer, intent(in) :: num
475  real (kind = 4) :: rieee(num)
476 
477  rieee(1:num) = transfer(cieee(1:num * 4), rieee, num)
478  call rdieee(rieee, a, num)
479 end subroutine rdieeec
480 
491 subroutine rdieee(rieee, a, num)
492  implicit none
493 
494  real(4), intent(in) :: rieee(num)
495  real, intent(out) :: a(num)
496  integer, intent(in) :: num
497 
498  integer(4) :: ieee
499  real, parameter :: two23 = scale(1.0, -23)
500  real, parameter :: two126 = scale(1.0, -126)
501  integer :: iexp, imant, isign, j
502  real :: sign, temp
503 
504  do j = 1, num
505  ! Transfer IEEE bit string to integer variable.
506  ieee = transfer(rieee(j), ieee)
507 
508  ! Extract sign bit, exponent, and mantissa.
509  isign = ibits(ieee, 31, 1)
510  iexp = ibits(ieee, 23, 8)
511  imant = ibits(ieee, 0, 23)
512  sign = 1.0
513  if (isign .eq. 1) sign = -1.0
514 
515  if (iexp .gt. 0 .and. iexp .lt. 255) then
516  temp = 2.0**(iexp - 127)
517  a(j) = sign * temp * (1.0 + (two23 * real(imant)))
518  elseif (iexp .eq. 0) then
519  if (imant .ne. 0) then
520  a(j) = sign * two126 * two23 * real(imant)
521  else
522  a(j) = sign * 0.0
523  endif
524  elseif (iexp .eq. 255) then
525  a(j) = sign * huge(a(j))
526  endif
527  enddo
528 end subroutine rdieee
529 
539 subroutine mkieee(a, rieee, num)
540  implicit none
541 
542  real(4), intent(in) :: a(num)
543  real(4), intent(out) :: rieee(num)
544  integer, intent(in) :: num
545 
546  integer(4) :: ieee
547  real, parameter :: two23 = scale(1.0,23)
548  real, parameter :: two126 = scale(1.0,126)
549  real :: alog2, atemp
550  integer :: iexp, imant, j, n
551 
552  alog2 = alog(2.0)
553 
554  do j = 1, num
555  ieee = 0
556  if (a(j) .eq. 0.) then
557  ieee = 0
558  rieee(j) = transfer(ieee, rieee(j))
559  cycle
560  endif
561 
562  ! Set Sign bit (bit 31 - leftmost bit).
563  if (a(j) .lt. 0.0) then
564  ieee = ibset(ieee, 31)
565  atemp = abs(a(j))
566  else
567  ieee = ibclr(ieee, 31)
568  atemp = a(j)
569  endif
570 
571  ! Determine exponent n with base 2.
572  if (atemp .ge. 1.0) then
573  n = 0
574  do while (2.0**(n+1) .le. atemp)
575  n = n + 1
576  enddo
577  else
578  n = -1
579  do while (2.0**n .gt. atemp )
580  n = n - 1
581  enddo
582  endif
583  iexp = n + 127
584  if (n .gt. 127) iexp = 255 ! overflow
585  if (n .lt. -127) iexp = 0
586  call mvbits(iexp, 0, 8, ieee, 23)
587 
588  ! Determine Mantissa.
589  if (iexp .ne. 255) then
590  if (iexp .ne. 0) then
591  atemp = (atemp / (2.0**n)) - 1.0
592  else
593  atemp = atemp * two126
594  endif
595  imant = nint(atemp * two23)
596  else
597  imant = 0
598  endif
599  ! set mantissa bits (bits 22-0).
600  call mvbits(imant, 0, 23, ieee, 0)
601 
602  ! Transfer IEEE bit string to real variable.
603  rieee(j) = transfer(ieee, rieee(j))
604  enddo
605 end subroutine mkieee
606 
subroutine g2_gbytesc(in, iout, iskip, nbits, nskip, n)
Extract arbitrary size big-endian integer values (up to 32 bits each) from a packed bit string.
Definition: g2bytes.F90:94
subroutine rdieee(rieee, a, num)
Copy array of 32-bit IEEE floating point values to local floating point representation.
Definition: g2bytes.F90:492
subroutine g2_sbytesc8(out, in, iskip, nbits, nskip, n)
Put arbitrary sized (up to 64 bits each) values into a packed bit string, taking the low order bits f...
Definition: g2bytes.F90:402
subroutine g2_gbytec(in, iout, iskip, nbits)
Extract one arbitrary size big-endian value (up to 32 bits) from a packed bit string into one element...
Definition: g2bytes.F90:21
subroutine mkieee(a, rieee, num)
Copy an array of real to an array of 32-bit IEEE floating points.
Definition: g2bytes.F90:540
subroutine g2_sbytescr(out, rin, iskip, nbits, nskip, n)
Put real values into a packed bit string in big-endian order.
Definition: g2bytes.F90:282
subroutine g2_gbytec8(in, iout, iskip, nbits)
Extract one arbitrary sized (up to 64-bits) values from a packed bit string, right justifying each va...
Definition: g2bytes.F90:152
subroutine g2_sbytec(out, in, iskip, nbits)
Put one arbitrary sized (up to 32 bits) value from an integer array, into a packed bit string,...
Definition: g2bytes.F90:238
subroutine g2_sbytec1(out, in, iskip, nbits)
Put one arbitrary sized (up to 32 bits) values from an integer scalar into a packed bit string,...
Definition: g2bytes.F90:260
subroutine rdieeec(cieee, a, num)
Copy array of 32-bit IEEE floating point values stored in char array to local floating point represen...
Definition: g2bytes.F90:470
subroutine g2_gbytesc8(in, iout, iskip, nbits, nskip, n)
Extract arbitrary sized (up to 64-bits) values from a packed bit string, right justifying each value ...
Definition: g2bytes.F90:173
subroutine g2_gbytescr(in, rout, iskip, nbits, nskip, n)
Extract big-endian floating-point values (32 bits each) from a packed bit string.
Definition: g2bytes.F90:68
subroutine g2_sbytec8(out, in, iskip, nbits)
Put one arbitrary sized (up to 64 bits) values into a packed bit string, taking the low order bits fr...
Definition: g2bytes.F90:380
subroutine g2_gbytec1(in, siout, iskip, nbits)
Extract one arbitrary size big-endian integer value (up to 32 bits) from a packed bit string into a s...
Definition: g2bytes.F90:43
subroutine g2_sbytesc(out, in, iskip, nbits, nskip, n)
Put arbitrary size (up to 32 bits each) integer values into a packed bit string in big-endian order.
Definition: g2bytes.F90:308