NCEPLIBS-g2  3.5.0
g2create.F90
Go to the documentation of this file.
1 
4 
53 subroutine gribcreate(cgrib, lcgrib, listsec0, listsec1, ierr)
54  implicit none
55 
56  character(len = 1), intent(inout) :: cgrib(lcgrib)
57  integer, intent(in) :: listsec0(*), listsec1(*)
58  integer, intent(in) :: lcgrib
59  integer, intent(out) :: ierr
60 
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
67 
68  interface
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
73  end subroutine g2_sbytec
74  subroutine g2_sbytec1(out, in, iskip, nbits)
75  character*1, intent(inout) :: out(*)
76  integer, intent(in) :: in
77  integer, intent(in) :: iskip, nbits
78  end subroutine g2_sbytec1
79  end interface
80 
81  ierr = 0
82 
83 #ifdef LOGGING
84  print *, 'gribcreate lcgrib ', lcgrib
85 #endif
86 
87  ! Currently handles only GRIB Edition 2.
88  if (listsec0(2) .ne. 2) then
89  print *, 'gribcreate: can only code GRIB edition 2.'
90  ierr = 1
91  return
92  endif
93 
94  ! Pack Section 0 - Indicator Section (except for total length of
95  ! GRIB message).
96  cgrib(1) = grib(1:1) ! Beginning of GRIB message
97  cgrib(2) = grib(2:2)
98  cgrib(3) = grib(3:3)
99  cgrib(4) = grib(4:4)
100  call g2_sbytec1(cgrib, zero, 32, 16) ! reserved for future use
101  call g2_sbytec(cgrib, listsec0(1), 48, 8) ! Discipline
102  call g2_sbytec(cgrib, listsec0(2), 56, 8) ! GRIB edition number
103  lensec0 = 16 ! bytes (octets)
104 
105  ! Pack Section 1 - Identification Section.
106  ibeg = lensec0 * 8 ! Calculate offset for beginning of section 1
107  iofst = ibeg + 32 ! leave space for length of section
108  call g2_sbytec1(cgrib, one, iofst, 8) ! Store section number ( 1 )
109  iofst = iofst + 8
110 
111  ! Pack up each input value in array listsec1 into the the
112  ! appropriate number of octets, which are specified in corresponding
113  ! entries in array mapsec1.
114  do i = 1, mapsec1len
115  nbits = mapsec1(i) * 8
116  call g2_sbytec(cgrib, listsec1(i), iofst, nbits)
117  iofst = iofst + nbits
118  enddo
119 
120  ! Calculate length of section 1 and store it in octets 1-4 of
121  ! section 1.
122  lensec1 = (iofst - ibeg) / 8
123  call g2_sbytec1(cgrib, lensec1, ibeg, 32)
124 
125  ! Put current byte total of message into Section 0.
126  call g2_sbytec1(cgrib, zero, 64, 32)
127  call g2_sbytec1(cgrib, lensec0 + lensec1, 96, 32)
128 end subroutine gribcreate
129 
196 subroutine addfield(cgrib, lcgrib, ipdsnum, ipdstmpl, ipdstmplen, &
197  coordlist, numcoord, idrsnum, idrstmpl, &
198  idrstmplen, fld, ngrdpts, ibmap, bmap, ierr)
199  use pdstemplates
200  use drstemplates
201  implicit none
202 
203  logical :: match
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)
214 
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)
233 
234  interface
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
239  end subroutine g2_gbytec
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
244  end subroutine g2_gbytec1
245  subroutine g2_gbytec81(in, siout, 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)
250  end subroutine g2_gbytec81
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
255  end subroutine g2_sbytec
256  subroutine g2_sbytec1(out, in, iskip, nbits)
257  character*1, intent(inout) :: out(*)
258  integer, intent(in) :: in
259  integer, intent(in) :: iskip, nbits
260  end subroutine g2_sbytec1
261  end interface
262 
263  allones = int(z'FFFFFFFF')
264  ierr = 0
265 
266  ! Check to see if beginning of GRIB message exists.
267  match = .true.
268  do i = 1, 4
269  if (cgrib(i) /= grib(i:i)) then
270  match = .false.
271  endif
272  enddo
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.'
276  ierr = 1
277  return
278  endif
279 
280  ! Get current length of GRIB message.
281  call g2_gbytec1(cgrib, lencurr, 96, 32)
282 
283  ! Check to see if GRIB message is already complete.
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.'
287  ierr = 2
288  return
289  endif
290 
291  ! Loop through all current sections of the GRIB message to find the
292  ! last section number.
293  issec3 = .false.
294  isprevbmap = .false.
295  len = 16 ! length of Section 0
296  do
297  ! Get number and length of next section
298  iofst = len * 8
299  call g2_gbytec1(cgrib, ilen, iofst, 32)
300  iofst = iofst + 32
301  call g2_gbytec1(cgrib, isecnum, iofst, 8)
302  iofst = iofst + 8
303 
304  ! Check if previous Section 3 exists and save location of
305  ! the section 3 in case needed later.
306  if (isecnum .eq. 3) then
307  issec3 = .true.
308  lpos3 = len + 1
309  lensec3 = ilen
310  endif
311 
312  ! Check if a previous defined bitmap exists
313  if (isecnum .eq. 6) then
314  call g2_gbytec1(cgrib, ibmprev, iofst, 8)
315  iofst = iofst + 8
316  if ((ibmprev .ge. 0) .and. (ibmprev .le. 253)) isprevbmap = .true.
317  endif
318  len = len + ilen
319 
320  ! Exit loop if last section reached
321  if (len .eq. lencurr) exit
322 
323  ! If byte count for each section does not match current
324  ! total length, then there is a problem.
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
329  ierr = 3
330  return
331  endif
332  enddo
333 
334  ! Sections 4 through 7 can only be added after section 3 or 7.
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.'
339  ierr = 4
340  return
341 
342  ! Sections 4 through 7 can only be added if section 3 was previously defined.
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.'
349  ierr = 6
350  return
351  endif
352 
353  ! Add Section 4 - Product Definition Section.
354  ibeg = lencurr * 8 ! Calculate offset for beginning of section 4
355  iofst = ibeg + 32 ! leave space for length of section
356  call g2_sbytec1(cgrib, four, iofst, 8) ! Store section number (4)
357  iofst = iofst + 8
358  call g2_sbytec1(cgrib, numcoord, iofst, 16) ! Store num of coordinate values
359  iofst = iofst + 16
360  call g2_sbytec1(cgrib, ipdsnum, iofst, 16) ! Store Prod Def Template num.
361  iofst = iofst + 16
362 
363  ! Get Product Definition Template.
364  call getpdstemplate(ipdsnum, mappdslen, mappds, needext, iret)
365  if (iret .ne. 0) then
366  ierr = 5
367  return
368  endif
369 
370  ! Extend the Product Definition Template, if necessary.
371  ! The number of values in a specific template may vary
372  ! depending on data specified in the "static" part of the
373  ! template.
374  if (needext) then
375  call extpdstemplate(ipdsnum, ipdstmpl, mappdslen, mappds)
376  endif
377 
378  ! Pack up each input value in array ipdstmpl into the
379  ! the appropriate number of octets, which are specified in
380  ! corresponding entries in array mappds.
381  do i = 1, mappdslen
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)
385  else
386  call g2_sbytec1(cgrib, one, iofst, 1)
387  call g2_sbytec1(cgrib, iabs(ipdstmpl(i)), iofst + 1, nbits - 1)
388  endif
389  iofst = iofst+nbits
390  enddo
391 
392  ! Add Optional list of vertical coordinate values
393  ! after the Product Definition Template, if necessary.
394  if (numcoord .ne. 0) then
395  do i = 1, numcoord
396  coordlist_4(i) = real(coordlist(i), 4)
397  end do
398  call mkieee(coordlist_4, coordieee, numcoord)
399  call g2_sbytescr(cgrib, coordieee, iofst, 32, 0, numcoord)
400  iofst = iofst + (32 * numcoord)
401  endif
402 
403  ! Calculate length of section 4 and store it in octets
404  ! 1-4 of section 4.
405  lensec4 = (iofst-ibeg)/8
406  call g2_sbytec1(cgrib, lensec4, ibeg, 32)
407 
408  ! Pack Data using appropriate algorithm
409 
410  ! Get Data Representation Template
411  call getdrstemplate(idrsnum, mapdrslen, mapdrs, needext, iret)
412  if (iret .ne. 0) then
413  ierr = 5
414  return
415  endif
416 
417  ! contract data field, removing data at invalid grid points,
418  ! if bit-map is provided with field.
419  if (ibmap .eq. 0 .OR. ibmap .eq. 254) then
420  allocate(pfld(max(2, ngrdpts)))
421  ndpts = 0;
422  do jj = 1, ngrdpts
423  intbmap(jj) = 0
424  if (bmap(jj)) then
425  intbmap(jj) = 1
426  ndpts = ndpts + 1
427  pfld(ndpts) = fld(jj);
428  endif
429  enddo
430  if (ndpts == 0 .and. ngrdpts > 0) then
431  pfld(1) = 0
432  endif
433  else
434  ndpts = ngrdpts;
435  pfld=>fld;
436  endif
437  lcpack = 0
438  nsize = ndpts*4
439  if (nsize .lt. minsize) nsize = minsize
440  allocate(cpack(nsize), stat = istat)
441  if (idrsnum.eq.0) then ! Simple Packing
442  call simpack(pfld, ndpts, idrstmpl, cpack, lcpack)
443  elseif (idrsnum.eq.2.or.idrsnum.eq.3) then ! Complex Packing
444  call cmplxpack(pfld, ndpts, idrsnum, idrstmpl, cpack, lcpack)
445  elseif (idrsnum.eq.50) then ! Sperical Harmonic Simple Packing
446  call simpack(pfld(2), ndpts-1, idrstmpl, cpack, lcpack)
447  tmpfld(1) = real(pfld(1), 4)
448  call mkieee(tmpfld, tmpre00, 1) ! ensure RE(0,0) value is IEEE format
449  re00 = tmpre00(1)
450  ire00 = transfer(re00, ire00)
451  idrstmpl(5) = ire00
452  elseif (idrsnum.eq.51) then ! Sperical Harmonic Complex Packing
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)
456  else
457  print *, 'addfield: Cannot pack DRT 5.51.'
458  ierr=9
459  return
460  endif
461 
462  elseif (idrsnum.eq.40 .OR. idrsnum.eq.40000) then ! JPEG2000 encoding
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
466  width=ndpts
467  height=1
468  elseif (width.eq.allones .OR. height.eq.allones) then
469  width=ndpts
470  height=1
471  elseif (ibits(iscan, 5, 1) .eq. 1) then ! Scanning mode: bit 3
472  itemp=width
473  width=height
474  height=itemp
475  endif
476  else
477  width=ndpts
478  height=1
479  endif
480  if (width<1 .or. height<1) then
481  ! Special case: bitmask off everywhere.
482  write(0, *) 'Warning: bitmask off everywhere.'
483  write(0, *) ' Pretend one point in jpcpack to avoid crash.'
484  width=1
485  height=1
486  endif
487  lcpack=nsize
488  !print *, 'w, h=', width, height
489  call jpcpack(pfld, width, height, idrstmpl, cpack, lcpack)
490 
491  elseif (idrsnum.eq.41 .OR. idrsnum.eq.40010) then ! PNG encoding
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
495  width=ndpts
496  height=1
497  elseif (width.eq.allones .OR. height.eq.allones) then
498  width=ndpts
499  height=1
500  elseif (ibits(iscan, 5, 1) .eq. 1) then ! Scanning mode: bit 3
501  itemp=width
502  width=height
503  height=itemp
504  endif
505  else
506  width=ndpts
507  height=1
508  endif
509  !print *, 'png size ', width, height
510  call pngpack(pfld, width, height, idrstmpl, cpack, lcpack)
511  !print *, 'png packed'
512  else
513  print *, 'addfield: Data Representation Template 5.', idrsnum, &
514  ' not yet implemented.'
515  ierr=7
516  return
517  endif
518  if (ibmap.eq.0 .OR. ibmap.eq.254) then
519  deallocate(pfld)
520  endif
521  if (lcpack .lt. 0) then
522  if (allocated(cpack))deallocate(cpack)
523  ierr=10
524  return
525  endif
526 
527  ! Add Section 5 - Data Representation Section
528  ibeg=iofst ! Calculate offset for beginning of section 5
529  iofst=ibeg + 32 ! leave space for length of section
530  call g2_sbytec1(cgrib, five, iofst, 8) ! Store section number (5)
531  iofst=iofst + 8
532  call g2_sbytec1(cgrib, ndpts, iofst, 32) ! Store num of actual data points
533  iofst=iofst + 32
534  call g2_sbytec1(cgrib, idrsnum, iofst, 16) ! Store Data Repr. Template num.
535  iofst = iofst + 16
536 
537  ! Pack up each input value in array idrstmpl into the
538  ! the appropriate number of octets, which are specified in
539  ! corresponding entries in array mapdrs.
540  do i=1, mapdrslen
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)
544  else
545  call g2_sbytec1(cgrib, one, iofst, 1)
546  call g2_sbytec1(cgrib, iabs(idrstmpl(i)), iofst+1, nbits-1)
547  endif
548  iofst = iofst + nbits
549  enddo
550 
551  ! Calculate length of section 5 and store it in octets
552  ! 1-4 of section 5.
553  lensec5 = (iofst - ibeg) / 8
554  call g2_sbytec1(cgrib, lensec5, ibeg, 32)
555 
556  ! Add Section 6 - Bit-Map Section.
557  ibeg = iofst ! Calculate offset for beginning of section 6
558  iofst = ibeg + 32 ! leave space for length of section
559  call g2_sbytec1(cgrib, six, iofst, 8) ! Store section number (6)
560  iofst=iofst + 8
561  call g2_sbytec1(cgrib, ibmap, iofst, 8) ! Store Bit Map indicator
562  iofst = iofst + 8
563 
564  ! Store bitmap, if supplied
565  if (ibmap.eq.0) then
566  call g2_sbytesc(cgrib, intbmap, iofst, 1, 0, ngrdpts) ! Store BitMap
567  iofst=iofst+ngrdpts
568  endif
569 
570  ! If specifying a previously defined bit-map, make sure
571  ! one already exists in the current GRIB message.
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.'
575  ierr=8
576  return
577  endif
578 
579  ! Calculate length of section 6 and store it in octets
580  ! 1-4 of section 6. Pad to end of octect, if necessary.
581  left = 8 - mod(iofst, 8)
582  if (left .ne. 8) then
583  call g2_sbytec1(cgrib, zero, iofst, left) ! Pad with zeros to fill Octet
584  iofst = iofst + left
585  endif
586  lensec6 = (iofst - ibeg) / 8
587  call g2_sbytec1(cgrib, lensec6, ibeg, 32)
588 
589  ! Add Section 7 - Data Section.
590  ibeg = iofst ! Calculate offset for beginning of section 7
591  iofst = ibeg + 32 ! leave space for length of section
592  call g2_sbytec1(cgrib, seven, iofst, 8) ! Store section number (7)
593  iofst = iofst + 8
594 
595  ! Store Packed Binary Data values, if non-constant field
596  if (lcpack .ne. 0) then
597  ioctet = iofst / 8
598  cgrib(ioctet + 1:ioctet + lcpack) = cpack(1:lcpack)
599  iofst = iofst + (8 * lcpack)
600  endif
601 
602  ! Calculate length of section 7 and store it in octets
603  ! 1-4 of section 7.
604  lensec7 = (iofst - ibeg) / 8
605  call g2_sbytec1(cgrib, lensec7, ibeg, 32)
606 
607  if (allocated(cpack)) deallocate(cpack)
608 
609  ! Update current byte total of message in Section 0.
610  newlen = lencurr + lensec4 + lensec5 + lensec6 + lensec7
611  call g2_sbytec1(cgrib, newlen, 96, 32)
612 
613  return
614 end subroutine addfield
615 
659 subroutine addgrid(cgrib, lcgrib, igds, igdstmpl, igdstmplen, &
660  ideflist, idefnum, ierr)
661  use gridtemplates
662  implicit none
663 
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
668 
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
674  logical needext
675  integer :: i, ilen, iret, isecnum, nbits
676 
677  interface
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
682  end subroutine g2_gbytec
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
687  end subroutine g2_gbytec1
688  subroutine g2_gbytec81(in, siout, 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)
693  end subroutine g2_gbytec81
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
698  end subroutine g2_sbytec
699  subroutine g2_sbytec1(out, in, iskip, nbits)
700  character*1, intent(inout) :: out(*)
701  integer, intent(in) :: in
702  integer, intent(in) :: iskip, nbits
703  end subroutine g2_sbytec1
704  end interface
705 
706  ierr = 0
707 
708 #ifdef LOGGING
709  print *, 'addgrid lcgrib ', lcgrib, ' igdstmplen ', igdstmplen, ' idefnum ', idefnum
710 #endif
711 
712  ! Check to see if beginning of GRIB message exists.
713  do i = 1, 4
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"')
719  print 10, cgrib(1:4)
720  ierr = 1
721  return
722  endif
723  enddo
724 
725  ! Get current length of GRIB message.
726  call g2_gbytec1(cgrib, lencurr, 96, 32)
727 
728  ! Check to see if GRIB message is already complete.
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', &
733  ' add new section.'
734  ierr = 2
735  return
736  endif
737 
738  ! Loop through all current sections of the GRIB message to find the
739  ! last section number.
740  len = 16 ! length of Section 0
741  do
742  ! Get length and section number of next section.
743  iofst = len * 8
744  call g2_gbytec1(cgrib, ilen, iofst, 32)
745  iofst = iofst + 32
746  call g2_gbytec1(cgrib, isecnum, iofst, 8)
747  len = len + ilen
748 
749  ! Exit loop if last section reached.
750  if (len .eq. lencurr) exit
751 
752  ! If byte count for each section doesn't match current
753  ! total length, then there is a problem.
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
758  ierr = 3
759  return
760  endif
761  enddo
762 
763  ! Section 3 can only be added after sections 1, 2 and 7.
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', &
767  ' 1, 2 or 7.'
768  print *, 'addgrid: Section ', isecnum, &
769  ' was the last found in given GRIB message.'
770  ierr = 4
771  return
772  endif
773 
774  ! Add Section 3 - Grid Definition Section.
775  ibeg = lencurr * 8 ! Calculate offset for beginning of section 3
776  iofst = ibeg + 32 ! leave space for length of section
777  call g2_sbytec1(cgrib, three, iofst, 8) ! Store section number (3)
778  iofst = iofst + 8
779  call g2_sbytec(cgrib, igds(1), iofst, 8) ! Store source of Grid def.
780  iofst = iofst + 8
781  call g2_sbytec(cgrib, igds(2), iofst, 32) ! Store number of data pts.
782  iofst = iofst + 32
783  call g2_sbytec(cgrib, igds(3), iofst, 8) ! Store number of extra octets.
784  iofst = iofst + 8
785  call g2_sbytec(cgrib, igds(4), iofst, 8) ! Store interp. of extra octets.
786  iofst = iofst + 8
787 
788  ! If Octet 6 is not equal to zero, Grid Definition Template may not
789  ! be supplied.
790  if (igds(1) .eq. 0) then
791  call g2_sbytec(cgrib, igds(5), iofst, 16) ! Store Grid Def Template num.
792  else
793  call g2_sbytec1(cgrib, 65535, iofst, 16) ! Store missing value as Grid Def Template num.
794  endif
795  iofst = iofst + 16
796 
797  ! Get Grid Definition Template.
798  if (igds(1) .eq. 0) then
799  call getgridtemplate(igds(5), mapgridlen, mapgrid, needext, &
800  iret)
801  if (iret .ne. 0) then
802  ierr = 5
803  return
804  endif
805 
806  ! Extend the Grid Definition Template, if necessary. The number
807  ! of values in a specific template may vary depending on data
808  ! specified in the "static" part of the template.
809  if (needext) then
810  call extgridtemplate(igds(5), igdstmpl, mapgridlen, &
811  mapgrid)
812  endif
813  else
814  mapgridlen = 0
815  endif
816 
817  ! Pack up each input value in array igdstmpl into the the
818  ! appropriate number of octets, which are specified in corresponding
819  ! entries in array mapgrid.
820  do i = 1, mapgridlen
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)
824  else
825  call g2_sbytec1(cgrib, one, iofst, 1)
826  call g2_sbytec1(cgrib, iabs(igdstmpl(i)), iofst + 1, nbits &
827  - 1)
828  endif
829  iofst = iofst + nbits
830  enddo
831 
832  ! If requested, insert optional list of numbers defining number of
833  ! points in each row or column. This is used for non regular grids.
834  if (igds(3) .ne. 0) then
835  nbits = igds(3) * 8
836  call g2_sbytesc(cgrib, ideflist, iofst, nbits, 0, idefnum)
837  iofst = iofst + (nbits * idefnum)
838  endif
839 
840  ! Calculate length of section 3 and store it in octets 1-4 of
841  ! section 3.
842  lensec3 = (iofst - ibeg) / 8
843  call g2_sbytec1(cgrib, lensec3, ibeg, 32)
844 
845  ! Update current byte total of message in Section 0.
846  call g2_sbytec1(cgrib, lencurr + lensec3, 96, 32)
847 end subroutine addgrid
848 
873 subroutine addlocal(cgrib, lcgrib, csec2, lcsec2, ierr)
874  implicit none
875 
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
880 
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
886 
887  ierr = 0
888 
889  ! Check to see if beginning of GRIB message exists.
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.'
895  ierr = 1
896  return
897  endif
898 
899  ! Get current length of GRIB message.
900  call g2_gbytec1(cgrib, lencurr, 96, 32)
901 
902  ! Check to see if GRIB message is already complete
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.'
906  ierr = 2
907  return
908  endif
909 
910  ! Loop through all current sections of the GRIB message to find the
911  ! last section number.
912  len = 16 ! length of Section 0
913  do
914  ! Get section number and length of next section.
915  iofst = len * 8
916  call g2_gbytec1(cgrib, ilen, iofst, 32)
917  iofst = iofst + 32
918  call g2_gbytec1(cgrib, isecnum, iofst, 8)
919  len = len + ilen
920  ! Exit loop if last section reached
921  if (len .eq. lencurr) exit
922 
923  ! If byte count for each section doesn't match current
924  ! total length, then there is a problem.
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
929  ierr = 3
930  return
931  endif
932  enddo
933 
934  ! Section 2 can only be added after sections 1 and 7.
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.'
938  ierr = 4
939  return
940  endif
941 
942  ! Add Section 2 - Local Use Section.
943  ibeg = lencurr * 8 ! Calculate offset for beginning of section 2
944  iofst = ibeg + 32 ! leave space for length of section
945  call g2_sbytec1(cgrib, two, iofst, 8) ! Store section number (2)
946  istart = lencurr + 5
947  cgrib(istart + 1:istart + lcsec2) = csec2(1:lcsec2)
948 
949  ! Calculate length of section 2 and store it in octets 1-4 of
950  ! section 2.
951  lensec2 = lcsec2 + 5 ! bytes
952  call g2_sbytec1(cgrib, lensec2, ibeg, 32)
953 
954  ! Update current byte total of message in Section 0.
955  call g2_sbytec1(cgrib, lencurr + lensec2, 96, 32)
956 
957 end subroutine addlocal
958 
979 subroutine gribend(cgrib, lcgrib, lengrib, ierr)
980  implicit none
981 
982  character(len = 1), intent(inout) :: cgrib(lcgrib)
983  integer, intent(in) :: lcgrib
984  integer, intent(out) :: lengrib, ierr
985 
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
991 
992  ierr = 0
993 
994  ! Check to see if beginning of GRIB message exists.
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.'
998  ierr = 1
999  return
1000  endif
1001 
1002  ! Get current length of GRIB message.
1003  call g2_gbytec1(cgrib, lencurr, 96, 32)
1004 
1005  ! Loop through all current sections of the GRIB message to
1006  ! find the last section number.
1007  len = 16 ! Length of Section 0
1008  do
1009  ! Get number and length of next section.
1010  iofst = len * 8
1011  call g2_gbytec1(cgrib, ilen, iofst, 32)
1012  iofst = iofst + 32
1013  call g2_gbytec1(cgrib, isecnum, iofst, 8)
1014  len = len + ilen
1015 
1016  ! Exit loop if last section reached.
1017  if (len .eq. lencurr) exit
1018 
1019  ! If byte count for each section doesn't match current total
1020  ! length, then there is a problem.
1021  if (len .gt. lencurr) then
1022  print *, 'gribend: Section byte counts don''t add ' &
1023  ,'to total.'
1024  print *, 'gribend: Sum of section byte counts = ', len
1025  print *, 'gribend: Total byte count in Section 0 = ', &
1026  lencurr
1027  ierr = 3
1028  return
1029  endif
1030  enddo
1031 
1032  ! Can only add End Section (Section 8) after Section 7.
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.'
1037  ierr = 4
1038  return
1039  endif
1040 
1041  ! Add Section 8 - End Section.
1042  cgrib(lencurr + 1:lencurr + 4) = c7777
1043 
1044  ! Update current byte total of message in Section 0.
1045  lengrib = lencurr + 4
1046  call g2_sbytec1(cgrib, lengrib, 96, 32)
1047 end subroutine gribend
subroutine cmplxpack(fld, ndpts, idrsnum, idrstmpl, cpack, lcpack)
Pack up a data field using a complex packing algorithm.
Definition: compack.F90:43
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:685
subroutine g2_sbytescr(out, rin, iskip, nbits, nskip, n)
Put real values into a packed bit string in big-endian order.
Definition: g2bytes.F90:368
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:306
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:337
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...
Definition: g2bytes.F90:209
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:52
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:402
subroutine addgrid(cgrib, lcgrib, igds, igdstmpl, igdstmplen, ideflist, idefnum, ierr)
Add a Grid Definition Section (Section 3) to a GRIB2 message.
Definition: g2create.F90:661
subroutine addlocal(cgrib, lcgrib, csec2, lcsec2, ierr)
Add a Local Use Section (Section 2) to a GRIB2 message.
Definition: g2create.F90:874
subroutine gribcreate(cgrib, lcgrib, listsec0, listsec1, ierr)
Initialize a new GRIB2 message and pack GRIB2 sections 0 (Indicator) and 1 (Identification).
Definition: g2create.F90:54
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.
Definition: g2create.F90:199
subroutine gribend(cgrib, lcgrib, lengrib, ierr)
Finalize a GRIB2 message after all grids and fields have been added.
Definition: g2create.F90:980
subroutine getdim(csec3, lcsec3, width, height, iscan)
Return the dimensions and scanning mode of a grid definition packed in the Grid Definition Section.
Definition: g2get.F90:1192
subroutine getpoly(csec3, lcsec3, jj, kk, mm)
Return the J, K, and M pentagonal resolution parameters specified in a GRIB2 Grid Definition Section ...
Definition: g2get.F90:1419
subroutine jpcpack(fld, width, height, idrstmpl, cpack, lcpack)
Pack a data field into a JPEG2000 code stream as defined in Data Representation Template 5....
Definition: g2jpc.F90:35
subroutine pngpack(fld, width, height, idrstmpl, cpack, lcpack)
Pack a data field into PNG image format, defined in [Data Representation Template 5....
Definition: g2png.F90:34
subroutine simpack(fld, ndpts, idrstmpl, cpack, lcpack)
Pack a data field using a simple packing algorithm.
Definition: g2sim.F90:27
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...
Definition: g2spec.F90:25
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.