NCEPLIBS-g2  3.5.0
g2gf.F90
Go to the documentation of this file.
1 
5 
34 subroutine putgb2(lugb, gfld, iret)
35  use grib_mod
36  implicit none
37 
38  integer, intent(in) :: lugb
39  type(gribfield), intent(in) :: gfld
40  integer, intent(out) :: iret
41 
42  character(len = 1), allocatable, dimension(:) :: cgrib
43  integer :: listsec0(2)
44  integer :: igds(5)
45  real :: coordlist
46  integer :: ilistopt
47  integer :: ierr, is, lcgrib, lengrib
48 
49  listsec0 = (/0, 2/)
50  igds = (/0, 0, 0, 0, 0/)
51  coordlist = 0.0
52  ilistopt = 0
53 
54  ! Figure out the maximum length of the GRIB2 message.
55  lcgrib = 16 + 21 + 4 ! Sections 0, 1, and 8.
56  ! Check for Section 2.
57  if (associated(gfld%local) .AND. gfld%locallen .gt. 0) then
58  lcgrib = lcgrib + gfld%locallen * 4
59  endif
60  ! Maximum size for Sections 3, 4, and 5 < 512 each.
61  lcgrib = lcgrib + 512 + 512 + 512
62  ! Is there a section 6?
63  if (gfld%ibmap .eq. 0) then
64  lcgrib = lcgrib + gfld%ngrdpts
65  endif
66  ! Section 7 holds the data.
67  lcgrib = lcgrib + gfld%ngrdpts * 4
68 
69 #ifdef LOGGING
70  print *, 'putgb2 lugb ', lugb, ' lcgrib ', lcgrib
71 #endif
72 
73  ! Allocate array for grib2 field.
74  allocate(cgrib(lcgrib), stat = is)
75  if (is .ne. 0) then
76  print *, 'putgb2: cannot allocate memory. ', is
77  iret = 2
78  endif
79 
80  ! Create new message.
81  listsec0(1) = gfld%discipline
82  listsec0(2) = gfld%version
83  if (associated(gfld%idsect)) then
84  call gribcreate(cgrib, lcgrib, listsec0, gfld%idsect, ierr)
85  if (ierr .ne. 0) then
86  write(6, *) 'putgb2: ERROR creating new GRIB2 field = ', ierr
87  endif
88  else
89  print *, 'putgb2: No Section 1 info available. '
90  iret = 10
91  deallocate(cgrib)
92  return
93  endif
94 
95  ! Add local use section to grib2 message.
96  if (associated(gfld%local) .AND. gfld%locallen .gt. 0) then
97  call addlocal(cgrib, lcgrib, gfld%local, gfld%locallen, ierr)
98  if (ierr .ne. 0) then
99  write(6, *) 'putgb2: ERROR adding local info = ', ierr
100  endif
101  endif
102 
103  ! Add grid to grib2 message.
104  igds(1) = gfld%griddef
105  igds(2) = gfld%ngrdpts
106  igds(3) = gfld%numoct_opt
107  igds(4) = gfld%interp_opt
108  igds(5) = gfld%igdtnum
109  if (associated(gfld%igdtmpl)) then
110  call addgrid(cgrib, lcgrib, igds, gfld%igdtmpl, gfld%igdtlen, &
111  ilistopt, gfld%num_opt, ierr)
112  if (ierr .ne. 0) then
113  write(6, *) 'putgb2: ERROR adding grid info = ', ierr
114  endif
115  else
116  print *, 'putgb2: No GDT info available. '
117  iret = 11
118  deallocate(cgrib)
119  return
120  endif
121 
122  ! Add data field to grib2 message.
123  if (associated(gfld%ipdtmpl) .AND. &
124  associated(gfld%idrtmpl) .AND. &
125  associated(gfld%fld)) then
126  call addfield(cgrib, lcgrib, gfld%ipdtnum, gfld%ipdtmpl, &
127  gfld%ipdtlen, coordlist, gfld%num_coord, &
128  gfld%idrtnum, gfld%idrtmpl, gfld%idrtlen, &
129  gfld%fld, gfld%ngrdpts, gfld%ibmap, gfld%bmap, &
130  ierr)
131  if (ierr .ne. 0) then
132  write(6, *) 'putgb2: ERROR adding data field = ', ierr
133  endif
134  else
135  print *, 'putgb2: Missing some field info. '
136  iret = 12
137  deallocate(cgrib)
138  return
139  endif
140 
141  ! Close grib2 message and write to file.
142  call gribend(cgrib, lcgrib, lengrib, ierr)
143  call wryte(lugb, lengrib, cgrib)
144 
145  deallocate(cgrib)
146 end subroutine putgb2
147 
210 subroutine gf_getfld(cgrib, lcgrib, ifldnum, unpack, expand, gfld, ierr)
211 
212  use grib_mod
213  implicit none
214 
215  character(len = 1), intent(in) :: cgrib(lcgrib)
216  integer, intent(in) :: lcgrib, ifldnum
217  logical, intent(in) :: unpack, expand
218  type(gribfield), intent(out) :: gfld
219  integer, intent(out) :: ierr
220 
221  character(len = 4), parameter :: grib = 'GRIB', c7777 = '7777'
222  character(len = 4) :: ctemp
223  real, pointer, dimension(:) :: newfld
224  integer:: listsec0(2), igds(5)
225  integer iofst, istart
226  logical*1, pointer, dimension(:) :: bmpsave
227  logical have3, have4, have5, have6, have7
228 
229  !implicit none additions
230  integer :: numfld, j, lengrib, ipos, lensec0, lensec
231  integer :: isecnum, jerr, n, numlocal
232 
233  interface
234  subroutine gf_unpack1(cgrib, lcgrib, iofst, ids, idslen, ierr)
235  character(len = 1), intent(in) :: cgrib(lcgrib)
236  integer, intent(in) :: lcgrib
237  integer, intent(inout) :: iofst
238  integer, pointer, dimension(:) :: ids
239  integer, intent(out) :: ierr, idslen
240  end subroutine gf_unpack1
241  subroutine gf_unpack2(cgrib, lcgrib, iofst, lencsec2, csec2, ierr)
242  character(len = 1), intent(in) :: cgrib(lcgrib)
243  integer, intent(in) :: lcgrib
244  integer, intent(inout) :: iofst
245  integer, intent(out) :: lencsec2
246  integer, intent(out) :: ierr
247  character(len = 1), pointer, dimension(:) :: csec2
248  end subroutine gf_unpack2
249  subroutine gf_unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, &
250  mapgridlen, ideflist, idefnum, ierr)
251  character(len = 1), intent(in) :: cgrib(lcgrib)
252  integer, intent(in) :: lcgrib
253  integer, intent(inout) :: iofst
254  integer, pointer, dimension(:) :: igdstmpl, ideflist
255  integer, intent(out) :: igds(5)
256  integer, intent(out) :: mapgridlen
257  integer, intent(out) :: ierr, idefnum
258  end subroutine gf_unpack3
259  subroutine gf_unpack4(cgrib, lcgrib, iofst, ipdsnum, ipdstmpl, &
260  mappdslen, coordlist, numcoord, ierr)
261  character(len = 1), intent(in) :: cgrib(lcgrib)
262  integer, intent(in) :: lcgrib
263  integer, intent(inout) :: iofst
264  real, pointer, dimension(:) :: coordlist
265  integer, pointer, dimension(:) :: ipdstmpl
266  integer, intent(out) :: ipdsnum
267  integer, intent(out) :: mappdslen
268  integer, intent(out) :: ierr, numcoord
269  end subroutine gf_unpack4
270  subroutine gf_unpack5(cgrib, lcgrib, iofst, ndpts, idrsnum, &
271  idrstmpl, mapdrslen, ierr)
272  character(len = 1), intent(in) :: cgrib(lcgrib)
273  integer, intent(in) :: lcgrib
274  integer, intent(inout) :: iofst
275  integer, intent(out) :: ndpts, idrsnum
276  integer, pointer, dimension(:) :: idrstmpl
277  integer, intent(out) :: mapdrslen
278  integer, intent(out) :: ierr
279  end subroutine gf_unpack5
280  subroutine gf_unpack6(cgrib, lcgrib, iofst, ngpts, ibmap, bmap, &
281  ierr)
282  character(len = 1), intent(in) :: cgrib(lcgrib)
283  integer, intent(in) :: lcgrib, ngpts
284  integer, intent(inout) :: iofst
285  integer, intent(out) :: ibmap
286  integer, intent(out) :: ierr
287  logical*1, pointer, dimension(:) :: bmap
288  end subroutine gf_unpack6
289  subroutine gf_unpack7(cgrib, lcgrib, iofst, igdsnum, igdstmpl, &
290  idrsnum, idrstmpl, ndpts, fld, ierr)
291  character(len = 1), intent(in) :: cgrib(lcgrib)
292  integer, intent(in) :: lcgrib, ndpts, idrsnum, igdsnum
293  integer, intent(inout) :: iofst
294  integer, pointer, dimension(:) :: idrstmpl, igdstmpl
295  integer, intent(out) :: ierr
296  real, pointer, dimension(:) :: fld
297  end subroutine gf_unpack7
298  end interface
299 
300  have3 = .false.
301  have4 = .false.
302  have5 = .false.
303  have6 = .false.
304  have7 = .false.
305  ierr = 0
306  numfld = 0
307  gfld%locallen = 0
308  nullify(gfld%list_opt, gfld%igdtmpl, gfld%ipdtmpl)
309  nullify(gfld%coord_list, gfld%idrtmpl, gfld%bmap, gfld%fld)
310 
311  ! Check for valid request number
312  if (ifldnum .le. 0) then
313  print *, 'gf_getfld: Request for field number ' &
314  ,'must be positive.'
315  ierr = 3
316  return
317  endif
318 
319  ! Check for beginning of GRIB message in the first 100 bytes
320  istart = 0
321  do j = 1, 100
322  ctemp = cgrib(j) // cgrib(j + 1) // cgrib(j + 2) // cgrib(j + 3)
323  if (ctemp .eq. grib) then
324  istart = j
325  exit
326  endif
327  enddo
328  if (istart .eq. 0) then
329  print *, 'gf_getfld: Beginning characters GRIB not found.'
330  ierr = 1
331  return
332  endif
333 
334  ! Unpack Section 0 - Indicator Section
335  iofst = 8 * (istart + 5)
336  call g2_gbytec(cgrib, listsec0(1), iofst, 8) ! Discipline
337  iofst = iofst + 8
338  call g2_gbytec(cgrib, listsec0(2), iofst, 8) ! GRIB edition number
339  iofst = iofst + 8
340  iofst = iofst + 32
341  call g2_gbytec(cgrib, lengrib, iofst, 32) ! Length of GRIB message
342  iofst = iofst + 32
343  lensec0 = 16
344  ipos = istart + lensec0
345 
346  ! Currently handles only GRIB Edition 2.
347  if (listsec0(2) .ne. 2) then
348  print *, 'gf_getfld: can only decode GRIB edition 2.'
349  ierr = 2
350  return
351  endif
352 
353  ! Loop through the remaining sections keeping track of the length of
354  ! each. Also keep the latest Grid Definition Section info. Unpack
355  ! the requested field number.
356  do
357  ! Check to see if we are at end of GRIB message
358  ctemp = cgrib(ipos) // cgrib(ipos + 1) // cgrib(ipos + 2) // &
359  cgrib(ipos + 3)
360  if (ctemp .eq. c7777) then
361  ipos = ipos + 4
362  ! If end of GRIB message not where expected, issue error
363  if (ipos.ne.(istart + lengrib)) then
364  print *, 'gf_getfld: "7777" found, but not ' &
365  ,'where expected.'
366  ierr = 4
367  return
368  endif
369  exit
370  endif
371 
372  ! Get length of Section and Section number
373  iofst = (ipos - 1) * 8
374  call g2_gbytec(cgrib, lensec, iofst, 32) ! Get Length of Section
375  iofst = iofst + 32
376  call g2_gbytec(cgrib, isecnum, iofst, 8) ! Get Section number
377  iofst = iofst + 8
378 
379  ! Check to see if section number is valid
380  if ((isecnum .lt. 1) .or. (isecnum .gt. 7)) then
381  print *, 'gf_getfld: Unrecognized Section Encountered = ', &
382  isecnum
383  ierr = 8
384  return
385  endif
386 
387  ! If found Section 1, decode elements in Identification Section.
388  if (isecnum .eq. 1) then
389  iofst = iofst - 40 ! reset offset to beginning of section
390  call gf_unpack1(cgrib, lcgrib, iofst, gfld%idsect, &
391  gfld%idsectlen, jerr)
392  if (jerr .ne. 0) then
393  ierr = 15
394  return
395  endif
396  endif
397 
398  ! If found Section 2, Grab local section. Save in case this is
399  ! the latest one before the requested field.
400  if (isecnum .eq. 2) then
401  iofst = iofst - 40 ! reset offset to beginning of section
402  if (associated(gfld%local)) deallocate(gfld%local)
403  call gf_unpack2(cgrib, lcgrib, iofst, gfld%locallen, &
404  gfld%local, jerr)
405  if (jerr .ne. 0) then
406  call gf_free(gfld)
407  ierr = 16
408  return
409  endif
410  endif
411 
412  ! If found Section 3, unpack the GDS info using the appropriate
413  ! template. Save in case this is the latest grid before the
414  ! requested field.
415  if (isecnum .eq. 3) then
416  iofst = iofst - 40 ! reset offset to beginning of section
417  if (associated(gfld%igdtmpl)) deallocate(gfld%igdtmpl)
418  if (associated(gfld%list_opt)) deallocate(gfld%list_opt)
419  call gf_unpack3(cgrib, lcgrib, iofst, igds, gfld%igdtmpl, &
420  gfld%igdtlen, gfld%list_opt, gfld%num_opt, jerr)
421  if (jerr .ne. 0) then
422  call gf_free(gfld)
423  ierr = 10
424  return
425  endif
426  have3 = .true.
427  gfld%griddef = igds(1)
428  gfld%ngrdpts = igds(2)
429  gfld%numoct_opt = igds(3)
430  gfld%interp_opt = igds(4)
431  gfld%igdtnum = igds(5)
432  endif
433 
434  ! If found Section 4, check to see if this field is the one
435  ! requested.
436  if (isecnum .eq. 4) then
437  numfld = numfld + 1
438  if (numfld .eq. ifldnum) then
439  gfld%discipline = listsec0(1)
440  gfld%version = listsec0(2)
441  gfld%ifldnum = ifldnum
442  gfld%unpacked = unpack
443  gfld%expanded = .false.
444  iofst = iofst-40 ! reset offset to beginning of section
445  call gf_unpack4(cgrib, lcgrib, iofst, gfld%ipdtnum, &
446  gfld%ipdtmpl, gfld%ipdtlen, gfld%coord_list, &
447  gfld%num_coord, jerr)
448  if (jerr .ne. 0) then
449  call gf_free(gfld)
450  ierr = 11
451  return
452  endif
453  have4 = .true.
454  endif
455  endif
456 
457  ! If found Section 5, check to see if this field is the one
458  ! requested.
459  if ((isecnum .eq. 5).and.(numfld .eq. ifldnum)) then
460  iofst = iofst-40 ! reset offset to beginning of section
461  call gf_unpack5(cgrib, lcgrib, iofst, gfld%ndpts, &
462  gfld%idrtnum, gfld%idrtmpl, gfld%idrtlen, jerr)
463  if (jerr .ne. 0) then
464  call gf_free(gfld)
465  ierr = 12
466  return
467  endif
468  have5 = .true.
469  endif
470 
471  ! If found Section 6, Unpack bitmap. Save in case this is the
472  ! latest bitmap before the requested field.
473  if (isecnum .eq. 6) then
474  if (unpack) then ! unpack bitmap
475  iofst = iofst - 40 ! reset offset to beginning of section
476  bmpsave => gfld%bmap ! save pointer to previous bitmap
477  call gf_unpack6(cgrib, lcgrib, iofst, gfld%ngrdpts, &
478  gfld%ibmap, gfld%bmap, jerr)
479  if (jerr .ne. 0) then
480  call gf_free(gfld)
481  ierr = 13
482  return
483  endif
484  have6 = .true.
485  if (gfld%ibmap .eq. 254) then ! use previously specified bitmap
486  if (associated(bmpsave)) then
487  gfld%bmap => bmpsave
488  else
489  print *, 'gf_getfld: Previous bit-map ' &
490  ,'specified, but none exists, '
491  call gf_free(gfld)
492  ierr = 17
493  return
494  endif
495  else ! get rid of it
496  if (associated(bmpsave)) deallocate(bmpsave)
497  endif
498  else ! do not unpack bitmap
499  call g2_gbytec(cgrib, gfld%ibmap, iofst, 8) ! Get BitMap Indicator
500  have6 = .true.
501  endif
502  endif
503 
504  ! If found Section 7, check to see if this field is the one
505  ! requested.
506  if ((isecnum .eq. 7) .and. (numfld .eq. ifldnum) .and. unpack) &
507  then
508  iofst = iofst - 40 ! reset offset to beginning of section
509  call gf_unpack7(cgrib, lcgrib, iofst, gfld%igdtnum, &
510  gfld%igdtmpl, gfld%idrtnum, &
511  gfld%idrtmpl, gfld%ndpts, &
512  gfld%fld, jerr)
513  if (jerr .ne. 0) then
514  call gf_free(gfld)
515  print *, 'gf_getfld: return from gf_unpack7 = ', jerr
516  ierr = 14
517  return
518  endif
519  have7 = .true.
520 
521  ! If bitmap is used with this field, expand data field
522  ! to grid, if possible.
523  if (gfld%ibmap .ne. 255 .AND. associated(gfld%bmap)) then
524  if (expand) then
525  allocate(newfld(gfld%ngrdpts))
526  n = 1
527  do j = 1, gfld%ngrdpts
528  if (gfld%bmap(j)) then
529  newfld(j) = gfld%fld(n)
530  n = n + 1
531  else
532  newfld(j) = 0.0
533  endif
534  enddo
535  deallocate(gfld%fld);
536  gfld%fld=>newfld;
537  gfld%expanded = .true.
538  else
539  gfld%expanded = .false.
540  endif
541  else
542  gfld%expanded = .true.
543  endif
544  endif
545 
546  ! Check to see if we read pass the end of the GRIB message and
547  ! missed the terminator string '7777'.
548  ipos = ipos + lensec ! Update beginning of section pointer
549  if (ipos .gt. (istart + lengrib)) then
550  print *, 'gf_getfld: "7777" not found at end ' &
551  ,'of GRIB message.'
552  call gf_free(gfld)
553  ierr = 7
554  return
555  endif
556  !
557  ! If unpacking requested, return when all sections have been
558  ! processed.
559  if (unpack .and. have3 .and. have4 .and. have5 .and. have6 &
560  .and. have7) return
561 
562  ! If unpacking is not requested, return when sections 3 through
563  ! 6 have been processed.
564  if ((.not. unpack) .and. have3 .and. have4 .and. have5 .and. &
565  have6) return
566  enddo
567 
568  ! If exited from above loop, the end of the GRIB message was reached
569  ! before the requested field was found.
570  print *, 'gf_getfld: GRIB message contained ', numlocal, &
571  ' different fields.'
572  print *, 'gf_getfld: The request was for the ', ifldnum, &
573  ' field.'
574  ierr = 6
575  call gf_free(gfld)
576 end subroutine gf_getfld
577 
584 subroutine gf_free(gfld)
585  use grib_mod
586  implicit none
587 
588  type(gribfield) :: gfld
589  integer :: is
590 
591  if (associated(gfld%idsect)) then
592  deallocate(gfld%idsect,stat=is)
593  endif
594  nullify(gfld%idsect)
595 
596  if (associated(gfld%local)) then
597  deallocate(gfld%local,stat=is)
598  endif
599  nullify(gfld%local)
600 
601  if (associated(gfld%list_opt)) then
602  deallocate(gfld%list_opt,stat=is)
603  endif
604  nullify(gfld%list_opt)
605 
606  if (associated(gfld%igdtmpl)) then
607  deallocate(gfld%igdtmpl,stat=is)
608  endif
609  nullify(gfld%igdtmpl)
610 
611  if (associated(gfld%ipdtmpl)) then
612  deallocate(gfld%ipdtmpl,stat=is)
613  endif
614  nullify(gfld%ipdtmpl)
615 
616  if (associated(gfld%coord_list)) then
617  deallocate(gfld%coord_list,stat=is)
618  endif
619  nullify(gfld%coord_list)
620 
621  if (associated(gfld%idrtmpl)) then
622  deallocate(gfld%idrtmpl,stat=is)
623  endif
624  nullify(gfld%idrtmpl)
625 
626  if (associated(gfld%bmap)) then
627  deallocate(gfld%bmap,stat=is)
628  endif
629  nullify(gfld%bmap)
630 
631  if (associated(gfld%fld)) then
632  deallocate(gfld%fld,stat=is)
633  endif
634  nullify(gfld%fld)
635 
636  return
637 end subroutine gf_free
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 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 gf_free(gfld)
Free memory that was used to store array values in derived type grib_mod::gribfield.
Definition: g2gf.F90:585
subroutine putgb2(lugb, gfld, iret)
Pack a field into a grib2 message and write that message to a file.
Definition: g2gf.F90:35
subroutine gf_getfld(cgrib, lcgrib, ifldnum, unpack, expand, gfld, ierr)
Return the Grid Definition, and Product Definition for a given data field.
Definition: g2gf.F90:211
subroutine gf_unpack1(cgrib, lcgrib, iofst, ids, idslen, ierr)
Unpack Section 1 (Identification Section) of a GRIB2 message, starting at octet 6 of that Section.
Definition: g2unpack.F90:42
subroutine gf_unpack2(cgrib, lcgrib, iofst, lencsec2, csec2, ierr)
Unpack Section 2 (Local Use Section) of a GRIB2 message.
Definition: g2unpack.F90:100
subroutine gf_unpack7(cgrib, lcgrib, iofst, igdsnum, igdstmpl, idrsnum, idrstmpl, ndpts, fld, ierr)
Unpack Section 7 (Data Section) of a GRIB2 message.
Definition: g2unpack.F90:650
subroutine gf_unpack5(cgrib, lcgrib, iofst, ndpts, idrsnum, idrstmpl, mapdrslen, ierr)
Unpack Section 5 (Data Representation Section) of a GRIB2 message, starting at octet 6 of that Sectio...
Definition: g2unpack.F90:467
subroutine gf_unpack4(cgrib, lcgrib, iofst, ipdsnum, ipdstmpl, mappdslen, coordlist, numcoord, ierr)
Unpack Section 4 (Product Definition Section) of a GRIB2 message, starting at octet 6 of that Section...
Definition: g2unpack.F90:334
subroutine gf_unpack3(cgrib, lcgrib, iofst, igds, igdstmpl, mapgridlen, ideflist, idefnum, ierr)
Unpack Section 3 (Grid Definition Section) of a GRIB2 message, starting at octet 6 of that Section.
Definition: g2unpack.F90:184
subroutine gf_unpack6(cgrib, lcgrib, iofst, ngpts, ibmap, bmap, ierr)
Unpack Section 6 (Bit-Map Section) of a GRIB2 message, starting at octet 6 of that Section.
Definition: g2unpack.F90:576
This Fortran module contains the declaration of derived type gribfield.
Definition: gribmod.F90:10