NCEPLIBS-g2  3.4.7
putgb2.F90
Go to the documentation of this file.
1 
6 
40 subroutine putgb2(lugb, gfld, iret)
41  use grib_mod
42  implicit none
43 
44  integer, intent(in) :: lugb
45  type(gribfield), intent(in) :: gfld
46  integer, intent(out) :: iret
47 
48  character(len = 1), allocatable, dimension(:) :: cgrib
49  integer :: listsec0(2)
50  integer :: igds(5)
51  real :: coordlist
52  integer :: ilistopt
53  integer :: ierr, is, lcgrib, lengrib
54 
55  listsec0 = (/0, 2/)
56  igds = (/0, 0, 0, 0, 0/)
57  coordlist = 0.0
58  ilistopt = 0
59 
60  ! Figure out the maximum length of the GRIB2 message.
61  lcgrib = 16 + 21 + 4 ! Sections 0, 1, and 8.
62  ! Check for Section 2.
63  if (associated(gfld%local) .AND. gfld%locallen .gt. 0) then
64  lcgrib = lcgrib + gfld%locallen * 4
65  endif
66  ! Maximum size for Sections 3, 4, and 5 < 512 each.
67  lcgrib = lcgrib + 512 + 512 + 512
68  ! Is there a section 6?
69  if (gfld%ibmap .eq. 0) then
70  lcgrib = lcgrib + gfld%ngrdpts
71  endif
72  ! Section 7 holds the data.
73  lcgrib = lcgrib + gfld%ngrdpts * 4
74 
75 #ifdef LOGGING
76  print *, 'putgb2 lugb ', lugb, ' lcgrib ', lcgrib
77 #endif
78 
79  ! Allocate array for grib2 field.
80  allocate(cgrib(lcgrib), stat = is)
81  if (is .ne. 0) then
82  print *, 'putgb2: cannot allocate memory. ', is
83  iret = 2
84  endif
85 
86  ! Create new message.
87  listsec0(1) = gfld%discipline
88  listsec0(2) = gfld%version
89  if (associated(gfld%idsect)) then
90  call gribcreate(cgrib, lcgrib, listsec0, gfld%idsect, ierr)
91  if (ierr .ne. 0) then
92  write(6, *) 'putgb2: ERROR creating new GRIB2 field = ', ierr
93  endif
94  else
95  print *, 'putgb2: No Section 1 info available. '
96  iret = 10
97  deallocate(cgrib)
98  return
99  endif
100 
101  ! Add local use section to grib2 message.
102  if (associated(gfld%local) .AND. gfld%locallen .gt. 0) then
103  call addlocal(cgrib, lcgrib, gfld%local, gfld%locallen, ierr)
104  if (ierr .ne. 0) then
105  write(6, *) 'putgb2: ERROR adding local info = ', ierr
106  endif
107  endif
108 
109  ! Add grid to grib2 message.
110  igds(1) = gfld%griddef
111  igds(2) = gfld%ngrdpts
112  igds(3) = gfld%numoct_opt
113  igds(4) = gfld%interp_opt
114  igds(5) = gfld%igdtnum
115  if (associated(gfld%igdtmpl)) then
116  call addgrid(cgrib, lcgrib, igds, gfld%igdtmpl, gfld%igdtlen, &
117  ilistopt, gfld%num_opt, ierr)
118  if (ierr .ne. 0) then
119  write(6, *) 'putgb2: ERROR adding grid info = ', ierr
120  endif
121  else
122  print *, 'putgb2: No GDT info available. '
123  iret = 11
124  deallocate(cgrib)
125  return
126  endif
127 
128  ! Add data field to grib2 message.
129  if (associated(gfld%ipdtmpl) .AND. &
130  associated(gfld%idrtmpl) .AND. &
131  associated(gfld%fld)) then
132  call addfield(cgrib, lcgrib, gfld%ipdtnum, gfld%ipdtmpl, &
133  gfld%ipdtlen, coordlist, gfld%num_coord, &
134  gfld%idrtnum, gfld%idrtmpl, gfld%idrtlen, &
135  gfld%fld, gfld%ngrdpts, gfld%ibmap, gfld%bmap, &
136  ierr)
137  if (ierr .ne. 0) then
138  write(6, *) 'putgb2: ERROR adding data field = ', ierr
139  endif
140  else
141  print *, 'putgb2: Missing some field info. '
142  iret = 12
143  deallocate(cgrib)
144  return
145  endif
146 
147  ! Close grib2 message and write to file.
148  call gribend(cgrib, lcgrib, lengrib, ierr)
149  call wryte(lugb, lengrib, cgrib)
150 
151  deallocate(cgrib)
152 end subroutine putgb2
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 given field and add them to a GRIB2 message.
Definition: addfield.f:75
subroutine addgrid(cgrib, lcgrib, igds, igdstmpl, igdstmplen, ideflist, idefnum, ierr)
Add a Grid Definition Section (Section 3) to a GRIB2 message.
Definition: addgrid.F90:51
subroutine addlocal(cgrib, lcgrib, csec2, lcsec2, ierr)
Add a Local Use Section (Section 2) to a GRIB2 message.
Definition: addlocal.F90:30
subroutine gribcreate(cgrib, lcgrib, listsec0, listsec1, ierr)
This subroutine initializes a new GRIB2 message and packs GRIB2 sections 0 (Indicator Section) and 1 ...
Definition: gribcreate.F90:55
subroutine gribend(cgrib, lcgrib, lengrib, ierr)
Finalize a GRIB2 message after all grids and fields have been added.
Definition: gribend.F90:27
This Fortran module contains the declaration of derived type gribfield.
Definition: gribmod.F90:10
subroutine putgb2(lugb, gfld, iret)
This subroutine packs a single field into a grib2 message and writes out that message to the file ass...
Definition: putgb2.F90:41