NCEPLIBS-g2  3.4.8
getgb2p.F90
Go to the documentation of this file.
1 
4 
92 subroutine getgb2p(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
93  extract, k, gribm, leng, iret)
94  use grib_mod
95  implicit none
96 
97  integer, intent(in) :: lugb, lugi, j, jdisc, jpdtn, jgdtn
98  integer, dimension(:) :: jids(*), jpdt(*), jgdt(*)
99  logical, intent(in) :: extract
100  integer, intent(out) :: k, iret, leng
101  character(len = 1), pointer, dimension(:) :: gribm
102 
103  type(gribfield) :: gfld
104  integer :: msk1, irgi, irgs, jk, lpos, msk2, mskp, nlen, nmess, nnum
105 
106  character(len = 1), pointer, dimension(:) :: cbuf
107  parameter(msk1 = 32000, msk2 = 4000)
108 
109  ! Declare interfaces (required for cbuf pointer).
110  interface
111  subroutine getg2i(lugi, cbuf, nlen, nnum, iret)
112  character(len = 1), pointer, dimension(:) :: cbuf
113  integer, intent(in) :: lugi
114  integer, intent(out) :: nlen, nnum, iret
115  end subroutine getg2i
116  subroutine getg2ir(lugb, msk1, msk2, mnum, cbuf, nlen, nnum, &
117  nmess, iret)
118  character(len = 1), pointer, dimension(:) :: cbuf
119  integer, intent(in) :: lugb, msk1, msk2, mnum
120  integer, intent(out) :: nlen, nnum, nmess, iret
121  end subroutine getg2ir
122  subroutine getgb2rp(lugb, cindex, extract, gribm, leng, iret)
123  integer, intent(in) :: lugb
124  character(len = 1), intent(in) :: cindex(*)
125  logical, intent(in) :: extract
126  integer, intent(out) :: leng, iret
127  character(len = 1), pointer, dimension(:) :: gribm
128  end subroutine getgb2rp
129  end interface
130 
131  ! Initialize the index information in cbuf.
132  irgi = 0
133  if (lugi .gt. 0) then
134  call getg2i(lugi, cbuf, nlen, nnum, irgi)
135  elseif (lugi .le. 0) then
136  mskp = 0
137  call getg2ir(lugb, msk1, msk2, mskp, cbuf, nlen, nnum, nmess, irgi)
138  endif
139  if (irgi .gt. 1) then
140  iret = 96
141  return
142  endif
143 
144  ! Find info from index and fill a grib_mod::gribfield variable.
145  call getgb2s(cbuf, nlen, nnum, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
146  jk, gfld, lpos, irgs)
147  if (irgs .ne. 0) then
148  iret = 99
149  call gf_free(gfld)
150  return
151  endif
152 
153  ! Extract grib message from file.
154  nullify(gribm)
155  call getgb2rp(lugb, cbuf(lpos:), extract, gribm, leng, iret)
156 
157  k = jk
158 
159  ! Free cbuf memory allocated in getg2i/getg2ir().
160  if (associated(cbuf)) deallocate(cbuf)
161 
162  call gf_free(gfld)
163 end subroutine getgb2p
subroutine getg2i(LUGI, CBUF, NLEN, NNUM, IRET)
Read a GRIB2 index file and return its contents.
Definition: getg2i.F90:59
subroutine getg2ir(LUGB, MSK1, MSK2, MNUM, CBUF, NLEN, NNUM, NMESS, IRET)
Generate an index record for a message in a GRIB2 file.
Definition: getg2ir.F90:38
subroutine getgb2p(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, extract, k, gribm, leng, iret)
Find and extract a GRIB2 message from a file.
Definition: getgb2p.F90:94
subroutine getgb2rp(lugb, cindex, extract, gribm, leng, iret)
Extract a grib message from a file given the index of the requested field.
Definition: getgb2rp.F90:39
subroutine getgb2s(CBUF, NLEN, NNUM, J, JDISC, JIDS, JPDTN, JPDT, JGDTN, JGDT, K, GFLD, LPOS, IRET)
Find information about a GRIB field from the index and fill a grib_mod::gribfield.
Definition: getgb2s.F90:80
subroutine gf_free(gfld)
Free memory that was used to store array values in derived type grib_mod::gribfield.
Definition: gf_free.F90:12
This Fortran module contains the declaration of derived type gribfield.
Definition: gribmod.F90:10