53 subroutine gribcreate(cgrib, lcgrib, listsec0, listsec1, ierr)
56 character(len = 1),
intent(inout) :: cgrib(lcgrib)
57 integer,
intent(in) :: listsec0(*), listsec1(*)
58 integer,
intent(in) :: lcgrib
59 integer,
intent(out) :: ierr
61 integer :: i, lensec1, nbits
62 character(len = 4),
parameter :: grib =
'GRIB'
63 integer,
parameter :: ZERO = 0, one = 1
64 integer,
parameter :: mapsec1len = 13
65 integer,
parameter :: mapsec1(mapsec1len) = (/ 2, 2, 1, 1, 1, 2, 1, 1, 1, 1, 1, 1, 1 /)
66 integer lensec0, iofst, ibeg
69 subroutine g2_sbytec(out, in, iskip, nbits)
70 character*1,
intent(inout) :: out(*)
71 integer,
intent(in) :: in(*)
72 integer,
intent(in) :: iskip, nbits
75 character*1,
intent(inout) :: out(*)
76 integer,
intent(in) :: in
77 integer,
intent(in) :: iskip, nbits
84 print *,
'gribcreate lcgrib ', lcgrib
88 if (listsec0(2) .ne. 2)
then
89 print *,
'gribcreate: can only code GRIB edition 2.'
101 call g2_sbytec(cgrib, listsec0(1), 48, 8)
102 call g2_sbytec(cgrib, listsec0(2), 56, 8)
115 nbits = mapsec1(i) * 8
116 call g2_sbytec(cgrib, listsec1(i), iofst, nbits)
117 iofst = iofst + nbits
122 lensec1 = (iofst - ibeg) / 8
127 call g2_sbytec1(cgrib, lensec0 + lensec1, 96, 32)
196 subroutine addfield(cgrib, lcgrib, ipdsnum, ipdstmpl, ipdstmplen, &
197 coordlist, numcoord, idrsnum, idrstmpl, &
198 idrstmplen, fld, ngrdpts, ibmap, bmap, ierr)
204 character(len=1),
intent(inout) :: cgrib(lcgrib)
205 integer,
intent(in) :: ipdsnum, ipdstmpl(*)
206 integer,
intent(in) :: idrsnum, numcoord, ipdstmplen, idrstmplen
207 integer,
intent(in) :: lcgrib, ngrdpts, ibmap
208 real,
intent(in) :: coordlist(numcoord)
209 real(kind = 4) :: coordlist_4(numcoord)
210 real,
target,
intent(in) :: fld(ngrdpts)
211 integer,
intent(out) :: ierr
212 integer,
intent(inout) :: idrstmpl(*)
213 logical*1,
intent(in) :: bmap(ngrdpts)
215 character(len=4),
parameter :: grib=
'GRIB', c7777=
'7777'
216 character(len=4):: ctemp
217 character(len=1),
allocatable :: cpack(:)
218 real,
pointer,
dimension(:) :: pfld
219 real(4) :: coordieee(numcoord), re00, tmpre00(1)
220 integer(4) :: ire00, allones
221 integer :: mappds(ipdstmplen), intbmap(ngrdpts), mapdrs(idrstmplen)
222 integer,
parameter :: zero=0, one=1, four=4, five=5, six=6, seven=7
223 integer,
parameter :: minsize=50000
224 integer iofst, ibeg, lencurr, len, mappdslen, mapdrslen, lpos3
225 integer width, height, ndpts
226 integer lensec3, lensec4, lensec5, lensec6, lensec7
227 logical issec3, needext, isprevbmap
228 integer :: nbits, newlen, nsize, lcpack, left
229 integer :: ibmprev, ilen, ioctet, iscan, isecnum, itemp
230 integer :: i, jj, kk, mm
231 integer :: iret, istat
232 real (kind = 4) :: tmpfld(1)
235 subroutine g2_gbytec(in, iout, iskip, nbits)
236 character*1,
intent(in) :: in(*)
237 integer,
intent(inout) :: iout(*)
238 integer,
intent(in) :: iskip, nbits
240 subroutine g2_gbytec1(in, siout, iskip, nbits)
241 character*1,
intent(in) :: in(*)
242 integer,
intent(inout) :: siout
243 integer,
intent(in) :: iskip, nbits
246 character*1,
intent(in) :: in(*)
247 integer (kind = 8),
intent(inout) :: siout
248 integer,
intent(in) :: iskip, nbits
249 integer (kind = 8) :: iout(1)
251 subroutine g2_sbytec(out, in, iskip, nbits)
252 character*1,
intent(inout) :: out(*)
253 integer,
intent(in) :: in(*)
254 integer,
intent(in) :: iskip, nbits
257 character*1,
intent(inout) :: out(*)
258 integer,
intent(in) :: in
259 integer,
intent(in) :: iskip, nbits
263 allones = int(z
'FFFFFFFF')
269 if (cgrib(i) /= grib(i:i))
then
273 if (.not. match)
then
274 print *,
'addfield: GRIB not found in given message.'
275 print *,
'addfield: Call to routine gribcreate required to initialize GRIB messge.'
284 ctemp = cgrib(lencurr-3) // cgrib(lencurr - 2) // cgrib(lencurr - 1) // cgrib(lencurr)
285 if (ctemp .eq. c7777)
then
286 print *,
'addfield: GRIB message already complete. Cannot add new section.'
306 if (isecnum .eq. 3)
then
313 if (isecnum .eq. 6)
then
316 if ((ibmprev .ge. 0) .and. (ibmprev .le. 253)) isprevbmap = .true.
321 if (len .eq. lencurr)
exit
325 if (len .gt. lencurr)
then
326 print *,
'addfield: Section byte counts don''t add to total.'
327 print *,
'addfield: Sum of section byte counts = ', len
328 print *,
'addfield: Total byte count in Section 0 = ', lencurr
335 if ((isecnum .ne. 3) .and. (isecnum .ne. 7))
then
336 print *,
'addfield: Sections 4-7 can only be added after Section 3 or 7.'
337 print *,
'addfield: Section ', isecnum,
' was the last found in', &
338 ' given GRIB message.'
343 elseif (.not. issec3)
then
344 print *,
'addfield: Sections 4-7 can only be added if Section', &
345 ' 3 was previously included.'
346 print *,
'addfield: Section 3 was not found in given GRIB message.'
347 print *,
'addfield: Call to routine addgrid required', &
348 ' to specify Grid definition.'
365 if (iret .ne. 0)
then
382 nbits = iabs(mappds(i)) * 8
383 if ((mappds(i) .ge. 0).or.(ipdstmpl(i) .ge. 0))
then
384 call g2_sbytec(cgrib, ipdstmpl(i), iofst, nbits)
387 call g2_sbytec1(cgrib, iabs(ipdstmpl(i)), iofst + 1, nbits - 1)
394 if (numcoord .ne. 0)
then
396 coordlist_4(i) = real(coordlist(i), 4)
398 call mkieee(coordlist_4, coordieee, numcoord)
399 call g2_sbytescr(cgrib, coordieee, iofst, 32, 0, numcoord)
400 iofst = iofst + (32 * numcoord)
405 lensec4 = (iofst-ibeg)/8
412 if (iret .ne. 0)
then
419 if (ibmap .eq. 0 .OR. ibmap .eq. 254)
then
420 allocate(pfld(max(2, ngrdpts)))
427 pfld(ndpts) = fld(jj);
430 if (ndpts == 0 .and. ngrdpts > 0)
then
439 if (nsize .lt. minsize) nsize = minsize
440 allocate(cpack(nsize), stat = istat)
441 if (idrsnum.eq.0)
then
442 call simpack(pfld, ndpts, idrstmpl, cpack, lcpack)
443 elseif (idrsnum.eq.2.or.idrsnum.eq.3)
then
444 call cmplxpack(pfld, ndpts, idrsnum, idrstmpl, cpack, lcpack)
445 elseif (idrsnum.eq.50)
then
446 call simpack(pfld(2), ndpts-1, idrstmpl, cpack, lcpack)
447 tmpfld(1) = real(pfld(1), 4)
448 call mkieee(tmpfld, tmpre00, 1)
450 ire00 = transfer(re00, ire00)
452 elseif (idrsnum.eq.51)
then
453 call getpoly(cgrib(lpos3), lensec3, jj, kk, mm)
454 if (jj.ne.0 .AND. kk.ne.0 .AND. mm.ne.0)
then
455 call specpack(pfld, ndpts, jj, kk, mm, idrstmpl, cpack, lcpack)
457 print *,
'addfield: Cannot pack DRT 5.51.'
462 elseif (idrsnum.eq.40 .OR. idrsnum.eq.40000)
then
463 if (ibmap.eq.255)
then
464 call getdim(cgrib(lpos3), lensec3, width, height, iscan)
465 if (width.eq.0 .OR. height.eq.0)
then
468 elseif (width.eq.allones .OR. height.eq.allones)
then
471 elseif (ibits(iscan, 5, 1) .eq. 1)
then
480 if (width<1 .or. height<1)
then
482 write(0, *)
'Warning: bitmask off everywhere.'
483 write(0, *)
' Pretend one point in jpcpack to avoid crash.'
489 call jpcpack(pfld, width, height, idrstmpl, cpack, lcpack)
491 elseif (idrsnum.eq.41 .OR. idrsnum.eq.40010)
then
492 if (ibmap.eq.255)
then
493 call getdim(cgrib(lpos3), lensec3, width, height, iscan)
494 if (width.eq.0 .OR. height.eq.0)
then
497 elseif (width.eq.allones .OR. height.eq.allones)
then
500 elseif (ibits(iscan, 5, 1) .eq. 1)
then
510 call pngpack(pfld, width, height, idrstmpl, cpack, lcpack)
513 print *,
'addfield: Data Representation Template 5.', idrsnum, &
514 ' not yet implemented.'
518 if (ibmap.eq.0 .OR. ibmap.eq.254)
then
521 if (lcpack .lt. 0)
then
522 if (
allocated(cpack))
deallocate(cpack)
541 nbits =iabs(mapdrs(i)) * 8
542 if ((mapdrs(i) .ge. 0) .or. (idrstmpl(i) .ge. 0))
then
543 call g2_sbytec(cgrib, idrstmpl(i), iofst, nbits)
546 call g2_sbytec1(cgrib, iabs(idrstmpl(i)), iofst+1, nbits-1)
548 iofst = iofst + nbits
553 lensec5 = (iofst - ibeg) / 8
566 call g2_sbytesc(cgrib, intbmap, iofst, 1, 0, ngrdpts)
572 if ((ibmap.eq.254).and.(.not.isprevbmap))
then
573 print *,
'addfield: Requested previously defined bitmap, ', &
574 ' but one does not exist in the current GRIB message.'
581 left = 8 - mod(iofst, 8)
582 if (left .ne. 8)
then
586 lensec6 = (iofst - ibeg) / 8
596 if (lcpack .ne. 0)
then
598 cgrib(ioctet + 1:ioctet + lcpack) = cpack(1:lcpack)
599 iofst = iofst + (8 * lcpack)
604 lensec7 = (iofst - ibeg) / 8
607 if (
allocated(cpack))
deallocate(cpack)
610 newlen = lencurr + lensec4 + lensec5 + lensec6 + lensec7
659 subroutine addgrid(cgrib, lcgrib, igds, igdstmpl, igdstmplen, &
660 ideflist, idefnum, ierr)
664 character(len = 1),
intent(inout) :: cgrib(lcgrib)
665 integer,
intent(in) :: igds(*), igdstmpl(*), ideflist(idefnum)
666 integer,
intent(in) :: lcgrib, idefnum, igdstmplen
667 integer,
intent(out) :: ierr
669 character(len = 4),
parameter :: grib =
'GRIB', c7777 =
'7777'
670 character(len = 4):: ctemp
671 integer:: mapgrid(igdstmplen)
672 integer,
parameter :: ONE = 1, three = 3
673 integer lensec3, iofst, ibeg, lencurr, len, mapgridlen
675 integer :: i, ilen, iret, isecnum, nbits
678 subroutine g2_gbytec(in, iout, iskip, nbits)
679 character*1,
intent(in) :: in(*)
680 integer,
intent(inout) :: iout(*)
681 integer,
intent(in) :: iskip, nbits
683 subroutine g2_gbytec1(in, siout, iskip, nbits)
684 character*1,
intent(in) :: in(*)
685 integer,
intent(inout) :: siout
686 integer,
intent(in) :: iskip, nbits
689 character*1,
intent(in) :: in(*)
690 integer (kind = 8),
intent(inout) :: siout
691 integer,
intent(in) :: iskip, nbits
692 integer (kind = 8) :: iout(1)
694 subroutine g2_sbytec(out, in, iskip, nbits)
695 character*1,
intent(inout) :: out(*)
696 integer,
intent(in) :: in(*)
697 integer,
intent(in) :: iskip, nbits
700 character*1,
intent(inout) :: out(*)
701 integer,
intent(in) :: in
702 integer,
intent(in) :: iskip, nbits
709 print *,
'addgrid lcgrib ', lcgrib,
' igdstmplen ', igdstmplen,
' idefnum ', idefnum
714 if (cgrib(i) /= grib(i:i))
then
715 print *,
'addgrid: GRIB not found in given message.'
716 print *,
'addgrid: Call to routine gribcreate required', &
717 ' to initialize GRIB messge.'
718 10
format(
'"', 4a1,
'" /= "GRIB"')
729 ctemp = cgrib(lencurr - 3) // cgrib(lencurr - 2) // cgrib(lencurr &
730 - 1) // cgrib(lencurr)
731 if (ctemp .eq. c7777)
then
732 print *,
'addgrid: GRIB message already complete. Cannot', &
750 if (len .eq. lencurr)
exit
754 if (len .gt. lencurr)
then
755 print *,
'addgrid: Section byte counts don''t add to total.'
756 print *,
'addgrid: Sum of section byte counts = ', len
757 print *,
'addgrid: Total byte count in Section 0 = ', lencurr
764 if ((isecnum .ne. 1) .and. (isecnum .ne. 2) .and. &
765 (isecnum .ne. 7))
then
766 print *,
'addgrid: Section 3 can only be added after Section', &
768 print *,
'addgrid: Section ', isecnum, &
769 ' was the last found in given GRIB message.'
781 call g2_sbytec(cgrib, igds(2), iofst, 32)
790 if (igds(1) .eq. 0)
then
791 call g2_sbytec(cgrib, igds(5), iofst, 16)
798 if (igds(1) .eq. 0)
then
801 if (iret .ne. 0)
then
821 nbits = iabs(mapgrid(i)) * 8
822 if ((mapgrid(i) .ge. 0) .or. (igdstmpl(i) .ge. 0))
then
823 call g2_sbytec(cgrib, igdstmpl(i), iofst, nbits)
826 call g2_sbytec1(cgrib, iabs(igdstmpl(i)), iofst + 1, nbits &
829 iofst = iofst + nbits
834 if (igds(3) .ne. 0)
then
836 call g2_sbytesc(cgrib, ideflist, iofst, nbits, 0, idefnum)
837 iofst = iofst + (nbits * idefnum)
842 lensec3 = (iofst - ibeg) / 8
846 call g2_sbytec1(cgrib, lencurr + lensec3, 96, 32)
873 subroutine addlocal(cgrib, lcgrib, csec2, lcsec2, ierr)
876 character(len = 1),
intent(inout) :: cgrib(lcgrib)
877 character(len = 1),
intent(in) :: csec2(lcsec2)
878 integer,
intent(in) :: lcgrib, lcsec2
879 integer,
intent(out) :: ierr
881 character(len = 4),
parameter :: grib =
'GRIB', c7777 =
'7777'
882 character(len = 4):: ctemp
883 integer,
parameter :: two = 2
884 integer :: lensec2, iofst, ibeg, lencurr, len
885 integer :: ilen, isecnum, istart
890 ctemp = cgrib(1) // cgrib(2) // cgrib(3) // cgrib(4)
891 if (ctemp .ne. grib)
then
892 print *,
'addlocal: GRIB not found in given message.'
893 print *,
'addlocal: Call to routine gribcreate required', &
894 ' to initialize GRIB messge.'
903 ctemp = cgrib(lencurr - 3) // cgrib(lencurr - 2) // cgrib(lencurr - 1) // cgrib(lencurr)
904 if (ctemp .eq. c7777)
then
905 print *,
'addlocal: GRIB message already complete. Cannot add new section.'
921 if (len .eq. lencurr)
exit
925 if (len .gt. lencurr)
then
926 print *,
'addlocal: Section byte counts don''t add to total.'
927 print *,
'addlocal: Sum of section byte counts = ', len
928 print *,
'addlocal: Total byte count in Section 0 = ', lencurr
935 if ((isecnum .ne. 1) .and. (isecnum .ne. 7))
then
936 print *,
'addlocal: Section 2 can only be added after Section 1 or Section 7.'
937 print *,
'addlocal: Section ', isecnum,
' was the last found in given GRIB message.'
947 cgrib(istart + 1:istart + lcsec2) = csec2(1:lcsec2)
955 call g2_sbytec1(cgrib, lencurr + lensec2, 96, 32)
979 subroutine gribend(cgrib, lcgrib, lengrib, ierr)
982 character(len = 1),
intent(inout) :: cgrib(lcgrib)
983 integer,
intent(in) :: lcgrib
984 integer,
intent(out) :: lengrib, ierr
986 integer ilen, isecnum
987 character(len = 4),
parameter :: grib =
'GRIB'
988 character(len = 1),
parameter :: c7777(4) = (/
'7',
'7',
'7',
'7' /)
989 character(len = 4):: ctemp
990 integer iofst, lencurr, len
995 ctemp = cgrib(1) // cgrib(2) // cgrib(3) // cgrib(4)
996 if (ctemp .ne. grib)
then
997 print *,
'gribend: GRIB not found in given message.'
1017 if (len .eq. lencurr)
exit
1021 if (len .gt. lencurr)
then
1022 print *,
'gribend: Section byte counts don''t add ' &
1024 print *,
'gribend: Sum of section byte counts = ', len
1025 print *,
'gribend: Total byte count in Section 0 = ', &
1033 if (isecnum .ne. 7)
then
1034 print *,
'gribend: Section 8 can only be added after Section 7.'
1035 print *,
'gribend: Section ', isecnum, &
1036 ' was the last found in',
' given GRIB message.'
1042 cgrib(lencurr + 1:lencurr + 4) = c7777
1045 lengrib = lencurr + 4
subroutine cmplxpack(fld, ndpts, idrsnum, idrstmpl, cpack, lcpack)
Pack up a data field using a complex packing algorithm.
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...
subroutine mkieee(a, rieee, num)
Copy an array of real to an array of 32-bit IEEE floating points.
subroutine g2_sbytescr(out, rin, iskip, nbits, nskip, n)
Put real values into a packed bit string in big-endian order.
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,...
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,...
subroutine g2_gbytec81(in, siout, iskip, nbits)
Extract one arbitrary size big-endian integer value (up to 64 bits) from a packed bit string into a s...
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...
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.
subroutine addgrid(cgrib, lcgrib, igds, igdstmpl, igdstmplen, ideflist, idefnum, ierr)
Add a Grid Definition Section (Section 3) to a GRIB2 message.
subroutine addlocal(cgrib, lcgrib, csec2, lcsec2, ierr)
Add a Local Use Section (Section 2) to a GRIB2 message.
subroutine gribcreate(cgrib, lcgrib, listsec0, listsec1, ierr)
Initialize a new GRIB2 message and pack GRIB2 sections 0 (Indicator) and 1 (Identification).
subroutine addfield(cgrib, lcgrib, ipdsnum, ipdstmpl, ipdstmplen, coordlist, numcoord, idrsnum, idrstmpl, idrstmplen, fld, ngrdpts, ibmap, bmap, ierr)
Pack up Sections 4 through 7 for a field and add them to a GRIB2 message.
subroutine gribend(cgrib, lcgrib, lengrib, ierr)
Finalize a GRIB2 message after all grids and fields have been added.
subroutine getdim(csec3, lcsec3, width, height, iscan)
Return the dimensions and scanning mode of a grid definition packed in the Grid Definition Section.
subroutine getpoly(csec3, lcsec3, jj, kk, mm)
Return the J, K, and M pentagonal resolution parameters specified in a GRIB2 Grid Definition Section ...
subroutine jpcpack(fld, width, height, idrstmpl, cpack, lcpack)
Pack a data field into a JPEG2000 code stream as defined in Data Representation Template 5....
subroutine pngpack(fld, width, height, idrstmpl, cpack, lcpack)
Pack a data field into PNG image format, defined in [Data Representation Template 5....
subroutine simpack(fld, ndpts, idrstmpl, cpack, lcpack)
Pack a data field using a simple packing algorithm.
subroutine specpack(fld, ndpts, JJ, KK, MM, idrstmpl, cpack, lcpack)
Pack a spectral data field using the complex packing algorithm for spherical harmonic data as defined...
Handles Data Representation Templates used in Section 5.
subroutine getdrstemplate(number, nummap, map, needext, iret)
Return DRS template information for a specified Data Representation Template.
This Fortran module contains info on all the available GRIB2 Grid Definition Templates used in [Secti...
subroutine getgridtemplate(number, nummap, map, needext, iret)
Get the grid template information for a specified Grid Definition Template.
subroutine extgridtemplate(number, list, nummap, map)
Generate the remaining octet map for a given Grid Definition Template, if required.
Information on all GRIB2 Product Definition Templates used in Section 4 - the Product Definition Sect...
subroutine extpdstemplate(number, list, nummap, map)
This subroutine generates the remaining octet map for a given Product Definition Template,...
subroutine getpdstemplate(number, nummap, map, needext, iret)
This subroutine returns PDS template information for a specified Product Definition Template.