NCEPLIBS-bufr
12.1.0
|
Decode character strings and integer values. More...
Go to the source code of this file.
Functions/Subroutines | |
recursive integer function | iupb (mbay, nbyt, nbit) |
Decode an integer value from within a specified number of bits of an integer array, starting with the first bit of a specified byte of the array. More... | |
recursive integer function | iupm (cbay, nbits) |
Decode an integer value from within a specified number of bits of a character string, starting with the first bit of the first byte of the string. More... | |
subroutine | up8 (nval, nbits, ibay, ibit) |
Decode an 8-byte integer value from within a specified number of bits of an integer array, starting at the bit immediately after a specified bit within the array. More... | |
subroutine | upb (nval, nbits, ibay, ibit) |
Decode an integer value from within a specified number of bits of an integer array, starting at the bit immediately after a specified bit within the array. More... | |
subroutine | upb8 (nval, nbits, ibit, ibay) |
Decode an 8-byte integer value from within a specified number of bits of an integer array, starting at the bit immediately after a specified bit within the array. More... | |
subroutine | upbb (nval, nbits, ibit, ibay) |
Decode an integer value from within a specified number of bits of an integer array, starting at the bit immediately after a specified bit within the array. More... | |
subroutine | upc (chr, nchr, ibay, ibit, cnvnull) |
Decode a character string from within a specified number of bytes of an integer array, starting at the bit immediately after a specified bit within the array. More... | |
real *8 function | ups (ival, node) |
Unpack a real*8 value from an integer by applying the proper scale and reference values. More... | |
recursive integer function iupb | ( | integer, dimension(*), intent(in) | mbay, |
integer, intent(in) | nbyt, | ||
integer, intent(in) | nbit | ||
) |
Decode an integer value from within a specified number of bits of an integer array, starting with the first bit of a specified byte of the array.
mbay | - Array containing encoded value |
nbyt | - Byte within mbay at whose first bit to begin decoding |
nbit | - Number of bits to be decoded, up to a maximum of 32 |
Definition at line 225 of file cidecode.F90.
Referenced by getlens(), iupbs01(), iupbs3(), rdmems(), rtrcptb(), stndrd(), upds3(), and wrdxtb().
recursive integer function iupm | ( | character*4, intent(in) | cbay, |
integer, intent(in) | nbits | ||
) |
Decode an integer value from within a specified number of bits of a character string, starting with the first bit of the first byte of the string.
This function is the logical inverse of subroutine ipkm().
cbay | - String |
nbits | - Number of bits from cbay to be decoded, up to a maximum of 32 |
Definition at line 264 of file cidecode.F90.
subroutine up8 | ( | integer*8, intent(out) | nval, |
integer, intent(in) | nbits, | ||
integer, dimension(*), intent(in) | ibay, | ||
integer, intent(inout) | ibit | ||
) |
Decode an 8-byte integer value from within a specified number of bits of an integer array, starting at the bit immediately after a specified bit within the array.
This subroutine is similar to subroutine upb8(), except that here ibit is both an input and an output argument, and the overall order of the arguments is different.
This subroutine is the logical inverse of subroutine pkb8().
ibay | - Array containing encoded value |
ibit | - Bit pointer within ibay
|
nbits | - Number of bits to be decoded |
nval | - Decoded value |
Definition at line 127 of file cidecode.F90.
References upb8().
subroutine upb | ( | integer, intent(out) | nval, |
integer, intent(in) | nbits, | ||
integer, dimension(*), intent(in) | ibay, | ||
integer, intent(inout) | ibit | ||
) |
Decode an integer value from within a specified number of bits of an integer array, starting at the bit immediately after a specified bit within the array.
This subroutine is similar to subroutine upbb(), except that here ibit is both an input and an output argument, and the overall order of the arguments is different.
ibay | - Array containing encoded value |
ibit | - Bit pointer within ibay
|
nbits | - Number of bits to be decoded |
nval | - Decoded value |
Definition at line 201 of file cidecode.F90.
References upbb().
Referenced by copysb(), iupb(), mvb(), rdcmps(), readsb(), stndrd(), ufbtab(), ufbtam(), upb8(), upc(), and writlc().
subroutine upb8 | ( | integer*8, intent(out) | nval, |
integer, intent(in) | nbits, | ||
integer, intent(in) | ibit, | ||
integer, dimension(*), intent(in) | ibay | ||
) |
Decode an 8-byte integer value from within a specified number of bits of an integer array, starting at the bit immediately after a specified bit within the array.
This subroutine is similar to subroutine up8(), except that here ibit is only an input argument, and the overall order of the arguments is different.
This subroutine will not work properly if nbits is greater than 64.
ibay | - Array containing encoded value |
ibit | - Bit within ibay after which to begin decoding nval |
nbits | - Number of bits to be decoded |
nval | - Decoded value |
Definition at line 79 of file cidecode.F90.
Referenced by rdtree(), ufbget(), ufbtab(), ufbtam(), and up8().
subroutine upbb | ( | integer, intent(out) | nval, |
integer, intent(in) | nbits, | ||
integer, intent(in) | ibit, | ||
integer, dimension(*), intent(in) | ibay | ||
) |
Decode an integer value from within a specified number of bits of an integer array, starting at the bit immediately after a specified bit within the array.
This subroutine is similar to subroutine upb(), except that here ibit is only an input argument, and the overall order of the arguments is different.
ibay | - Array containing encoded value |
ibit | - Bit within ibay after which to begin decoding nval |
nbits | - Number of bits to be decoded |
nval | - Decoded value |
Definition at line 153 of file cidecode.F90.
subroutine upc | ( | character*(*), intent(out) | chr, |
integer, intent(in) | nchr, | ||
integer, dimension(*), intent(in) | ibay, | ||
integer, intent(inout) | ibit, | ||
logical, intent(in) | cnvnull | ||
) |
Decode a character string from within a specified number of bytes of an integer array, starting at the bit immediately after a specified bit within the array.
chr | - Decoded string |
nchr | - Number of bytes of ibay from within which to decode chr (i.e. the number of characters in chr) |
ibay | - Array from which to decode chr |
ibit | - Bit pointer within ibay
|
cnvnull | - .true. if null characters in ibay should be converted to blanks within chr; .false. otherwise |
Definition at line 25 of file cidecode.F90.
References upb().
Referenced by rdcmps(), rdtree(), readlc(), stbfdx(), stndrd(), ufbget(), ufbtab(), ufbtam(), and wrcmps().
real*8 function ups | ( | integer*8, intent(in) | ival, |
integer, intent(in) | node | ||
) |
Unpack a real*8 value from an integer by applying the proper scale and reference values.
Normally the scale and reference values are obtained from index node of the internal jump/link table arrays isc(*) and irf(*); however, the reference value in irf(*) will be overridden if a 2-03 operator is in effect for this node.
This function is the logical inverse of function ipks().
ival | - Packed BUFR integer |
node | - Index into internal jump/link tables |
Definition at line 318 of file cidecode.F90.
References moda_tables::ibt, moda_nrv203::ienrv, moda_nrv203::inodnrv, moda_tables::irf, moda_tables::isc, moda_nrv203::isnrv, moda_nrv203::nnrv, moda_nrv203::nrv, moda_tables::tag, and moda_nrv203::tagnrv.
Referenced by rdcmps(), rdtree(), ufbget(), ufbtab(), and ufbtam().