NCEPLIBS-g2  3.4.7
getidx.F90
Go to the documentation of this file.
1 
5 
43 subroutine getidx(lugb, lugi, cindex, nlen, nnum, iret)
44  implicit none
45 
46  integer, intent(in) :: lugb, lugi
47  integer, intent(out) :: nlen, nnum, iret
48  character(len = 1), pointer, dimension(:) :: cindex
49  integer, parameter :: maxidx = 10000
50  integer, parameter :: msk1 = 32000, msk2 = 4000
51 
52  integer :: lux
53  integer :: irgi, mskp, nmess, i
54 
55  type gindex
56  integer :: nlen
57  integer :: nnum
58  character(len = 1), pointer, dimension(:) :: cbuf
59  end type gindex
60 
61  type(gindex), save :: idxlist(10000)
62 
63  data lux/0/
64 
65  ! declare interfaces (required for cbuf pointer)
66  interface
67  subroutine getg2i(lugi, cbuf, nlen, nnum, iret)
68  character(len = 1), pointer, dimension(:) :: cbuf
69  integer, intent(in) :: lugi
70  integer, intent(out) :: nlen, nnum, iret
71  end subroutine getg2i
72  subroutine getg2ir(lugb, msk1, msk2, mnum, cbuf, nlen, nnum, &
73  nmess, iret)
74  character(len = 1), pointer, dimension(:) :: cbuf
75  integer, intent(in) :: lugb, msk1, msk2, mnum
76  integer, intent(out) :: nlen, nnum, nmess, iret
77  end subroutine getg2ir
78  end interface
79 
80  ! Free all associated memory and exit.
81  if (lugb .eq. 0) then
82  !print *, 'getidx: Freeing all memory'
83  do i = 1, 10000
84  if (associated(idxlist(i)%cbuf)) then
85  !print *, 'deallocating ', loc(idxlist(i)%cbuf)
86  deallocate(idxlist(i)%cbuf)
87  nullify(idxlist(i)%cbuf)
88  endif
89  end do
90  iret = 0
91  return
92  endif
93 
94  ! determine whether index buffer needs to be initialized
95  lux = 0
96  iret = 0
97  if (lugb .le. 0 .or. lugb .gt. 9999) then
98  print *, ' file unit number out of range'
99  print *, ' use unit numbers in range: 0 - 9999 '
100  iret = 90
101  return
102  endif
103  if (lugi .eq. lugb) then ! force regeneration of index from grib2 file
104  if (associated(idxlist(lugb)%cbuf)) &
105  deallocate(idxlist(lugb)%cbuf)
106  !print *, 'Force regeneration'
107  nullify(idxlist(lugb)%cbuf)
108  idxlist(lugb)%nlen = 0
109  idxlist(lugb)%nnum = 0
110  lux = 0
111  endif
112 
113  if (lugi .lt. 0) then ! force re-read of index from indexfile
114  ! associated with unit abs(lugi)
115  if (associated(idxlist(lugb)%cbuf)) &
116  deallocate(idxlist(lugb)%cbuf)
117  !print *, 'Force re-read'
118  nullify(idxlist(lugb)%cbuf)
119  idxlist(lugb)%nlen = 0
120  idxlist(lugb)%nnum = 0
121  lux = abs(lugi)
122  endif
123 
124  ! check if index already exists in memory
125  if (associated(idxlist(lugb)%cbuf)) then
126  !print *, 'Index exists in memory!'
127  cindex => idxlist(lugb)%cbuf
128  nlen = idxlist(lugb)%nlen
129  nnum = idxlist(lugb)%nnum
130  return
131  endif
132 
133  irgi = 0
134  if (lux .gt. 0) then
135  call getg2i(lux, idxlist(lugb)%cbuf, nlen, nnum, irgi)
136  elseif (lux .le. 0) then
137  mskp = 0
138  call getg2ir(lugb, msk1, msk2, mskp, idxlist(lugb)%cbuf, &
139  nlen, nnum, nmess, irgi)
140  endif
141  if (irgi .eq. 0) then
142  cindex => idxlist(lugb)%cbuf
143  idxlist(lugb)%nlen = nlen
144  idxlist(lugb)%nnum = nnum
145  else
146  nlen = 0
147  nnum = 0
148  print *, ' error reading index file '
149  iret = 96
150  return
151  endif
152 end subroutine getidx
153 
160 subroutine gf_finalize(iret)
161  implicit none
162 
163  integer, intent(out) :: iret
164  character(len = 1), pointer, dimension(:) :: cindex
165  integer :: nlen, nnum
166 
167  ! Declare interfaces (required for cbuf pointer).
168  interface
169  subroutine getidx(lugb,lugi,cbuf,nlen,nnum,irgi)
170  character(len=1),pointer,dimension(:) :: cbuf
171  integer,intent(in) :: lugb,lugi
172  integer,intent(out) :: nlen,nnum,irgi
173  end subroutine getidx
174  end interface
175 
176  ! Call getidx with 0 for the first parameter, ensuring that the
177  ! internal memory is freed.
178  call getidx(0, 0, cindex, nlen, nnum, iret)
179 
180 end subroutine gf_finalize
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 getidx(lugb, lugi, cindex, nlen, nnum, iret)
Find, read or generate a GRIB2 index for the GRIB2 file associated with unit lugb.
Definition: getidx.F90:44
subroutine gf_finalize(iret)
Free all memory associated with the library.
Definition: getidx.F90:161