NCEPLIBS-g2 4.0.0
Loading...
Searching...
No Matches
cnv22.F90
Go to the documentation of this file.
1
5
37subroutine cnv22(ifl1,ifl2,ipack,usemiss,imiss,table_ver)
38
39 use grib_mod
40 use params
41 use re_alloc
42
43 integer,intent(in) :: ifl1,ifl2,ipack
44 logical,intent(in) :: usemiss
45
46 CHARACTER(len=1),pointer,dimension(:) :: cgrib
47 CHARACTER(len=8) :: ctemp
48 type(gribfield) :: gfld,prevfld
49 integer,dimension(200) :: jids,jpdt,jgdt
50 integer :: listsec0(2)=(/0,2/)
51 integer :: igds(5)=(/0,0,0,0,0/),previgds(5)
52 integer :: idrstmpl(200)
53 integer :: currlen=1000000
54 integer :: table_ver
55 logical :: unpack=.true.
56 logical :: open_grb=.false.
57 !
58 ! --- Initialize Variables ---
59 !
60 gfld%idsect => null()
61 gfld%local => null()
62 gfld%list_opt => null()
63 gfld%igdtmpl => null()
64 gfld%ipdtmpl => null()
65 gfld%coord_list => null()
66 gfld%idrtmpl => null()
67 gfld%bmap => null()
68 gfld%fld => null()
69
70 allocate(cgrib(currlen))
71 ifli1=0
72 jdisc=-1
73 jids=-9999
74 jpdt=-9999
75 jgdt=-9999
76 jpdtn=-1
77 jgdtn=-1
78 !
79 npoints=0
80 icount=0
81 jskp=0
82 do
83 prevfld=gfld
84 call getgb2(ifl1,ifli1,jskp,jdisc,jids,jpdtn,jpdt,jgdtn,jgdt, &
85 unpack,jskp,gfld,iret)
86 if (iret.ne.0) then
87 if (iret.eq.99) exit
88 print *,' getgb2 error = ',iret
89 cycle
90 !call errexit(17)
91 endif
92 icount=icount+1
93 !
94 ! Ensure that cgrib array is large enough
95 !
96 if (gfld%ifldnum == 1) then ! start new GRIB2 message
97 npoints=gfld%ngrdpts
98 else
99 npoints=npoints+gfld%ngrdpts
100 endif
101 newlen=npoints*4
102 if (newlen.gt.currlen) then
103 !if (allocated(cgrib)) deallocate(cgrib)
104 !allocate(cgrib(newlen),stat=is)
105 call realloc(cgrib,currlen,newlen,is)
106 currlen=newlen
107 endif
108 !
109 ! Start new GRIB2 message, if necessary.
110 ! May have to finish the current message though.
111 !
112 if (gfld%ifldnum == 1) then ! start new GRIB2 message
113 if (open_grb) then ! close previous GRIB2 message first
114 call gribend(cgrib,lcgrib,lengrib,ierr)
115 if (ierr.ne.0) then
116 write(6,*) ' ERROR ending new GRIB2 message = ',ierr
117 cycle
118 endif
119 open_grb=.false.
120 call wryte(ifl2,lengrib,cgrib)
121 endif
122 !
123 ! Create new GRIB Message
124 !
125 listsec0(1)=gfld%discipline
126 listsec0(2)=gfld%version
127 gfld%idsect(3) = table_ver
128 call gribcreate(cgrib,lcgrib,listsec0,gfld%idsect,ierr)
129 if (ierr.ne.0) then
130 write(6,*) ' ERROR creating new GRIB2 field = ',ierr
131 cycle
132 endif
133 open_grb=.true.
134 endif
135 !
136 ! Add grid to GRIB message, if previous grid in same
137 ! message is not the same.
138 !
139 previgds=igds
140 igds(1)=gfld%griddef
141 igds(2)=gfld%ngrdpts
142 igds(3)=gfld%numoct_opt
143 igds(4)=gfld%interp_opt
144 igds(5)=gfld%igdtnum
145 if (.NOT. associated(gfld%list_opt)) &
146 allocate(gfld%list_opt(1))
147 if (gfld%ifldnum == 1) then ! add grid to GRIB2 message
148 call addgrid(cgrib,lcgrib,igds,gfld%igdtmpl,gfld%igdtlen, &
149 gfld%list_opt,gfld%num_opt,ierr)
150 else ! check if previous grid is the same as the current
151 if (gfld%igdtlen.ne.prevfld%igdtlen .OR. &
152 gfld%num_opt.ne.prevfld%num_opt .OR. &
153 any(igds.ne.previgds) .OR. &
154 any(gfld%igdtmpl(1:gfld%igdtlen).NE. &
155 prevfld%igdtmpl(1:prevfld%igdtlen)) .OR. &
156 any(gfld%list_opt(1:gfld%num_opt).NE. &
157 prevfld%list_opt(1:prevfld%num_opt))) then
158 call addgrid(cgrib,lcgrib,igds,gfld%igdtmpl,gfld%igdtlen, &
159 gfld%list_opt,gfld%num_opt,ierr)
160 endif
161 endif
162 if (ierr.ne.0) then
163 write(6,*) ' ERROR adding GRIB2 grid = ',ierr
164 cycle
165 endif
166 call gf_free(prevfld)
167 idrstmpl=0
168 !
169 ! if usemiss is specified, change any bitmaps to
170 ! missing value management for DRTs 5.2 and 5.3.
171 ! OR carry on missing value management for fields
172 ! already using it.
173 !
174 if (usemiss .AND. &
175 (ipack.eq.2 .OR. ipack.eq.31 .OR. ipack.eq.32)) then
176 if (gfld%ibmap.eq.0 .OR. gfld%ibmap.eq.254) then
177 ! change bit-map to missing value mngmt.
178 gfld%ibmap=255
179 rmiss=minval(gfld%fld(1:gfld%ngrdpts))
180 if (rmiss .lt. -9999.0) then
181 rmiss=rmiss*10.0
182 else
183 rmiss=-9999.0
184 endif
185 do i=1,gfld%ngrdpts
186 if (.NOT. gfld%bmap(i)) then
187 gfld%fld(i)=rmiss
188 gfld%bmap(i)=.true.
189 endif
190 enddo
191 idrstmpl(7)=imiss ! Primary missing values
192 call mkieee(rmiss,idrstmpl(8),1)
193 elseif (gfld%idrtnum.EQ.2 .OR. gfld%idrtnum.EQ.3) then
194 idrstmpl(7)=gfld%idrtmpl(7) ! Missing value mgmt
195 idrstmpl(8)=gfld%idrtmpl(8) ! Primary missing value
196 idrstmpl(9)=gfld%idrtmpl(9) ! Secondary missing value
197 endif
198 endif
199 !
200 ! If converting from a field using missing value management
201 ! in DRTs 5.2 and 5.3 to a DRT that does not support missing
202 ! values, convert missings to a bitmap.
203 !
204 if ((.NOT. usemiss) .AND. &
205 (gfld%idrtnum.EQ.2 .OR. gfld%idrtnum.EQ.3) .AND. &
206 (gfld%idrtmpl(7).EQ.1 .OR. gfld%idrtmpl(7).EQ.2)) then
207 call rdieee(gfld%idrtmpl(8),rmissp,1)
208 if (gfld%idrtmpl(7) .EQ. 2) then
209 call rdieee(gfld%idrtmpl(9),rmisss,1)
210 else
211 rmisss=rmissp
212 endif
213 allocate(gfld%bmap(gfld%ngrdpts))
214 do j=1,gfld%ngrdpts
215 if (gfld%fld(j).EQ.rmissp .OR. &
216 gfld%fld(j).EQ.rmisss) then
217 gfld%bmap(j)=.false.
218 else
219 gfld%bmap(j)=.true.
220 endif
221 enddo
222 gfld%ibmap=0
223 idrstmpl(7)=0
224 idrstmpl(8)=0
225 idrstmpl(9)=0
226 endif
227 !
228 ! Add field to GRIB message
229 !
230 ! Set DRT info (packing info)
231 if (ipack.eq.0) then
232 idrsnum=0
233 elseif (ipack.eq.2) then
234 idrsnum=2
235 idrstmpl(6)=1
236 elseif (ipack.eq.31.OR.ipack.eq.32) then
237 idrsnum=ipack/10
238 idrstmpl(6)=1
239 idrstmpl(17)=mod(ipack,10) ! order of s.d.
240 elseif (ipack.eq.40 .OR. ipack.eq.41 .OR. &
241 ipack.eq.40000 .OR. ipack.eq.40010) then
242 idrsnum=ipack
243 idrstmpl(6)=0
244 idrstmpl(7)=255
245 else
246 idrsnum=3
247 idrstmpl(17)=1 ! order of s.d.
248 idrstmpl(6)=1 ! general group split
249 ctemp=param_get_abbrev(gfld%discipline,gfld%ipdtmpl(1), &
250 gfld%ipdtmpl(2))
251 if (ctemp.eq.'A PCP ') idrsnum=2
252 endif
253 idrstmpl(2)=gfld%idrtmpl(2)
254 idrstmpl(3)=gfld%idrtmpl(3)
255 if (.NOT. associated(gfld%coord_list)) &
256 allocate(gfld%coord_list(1))
257 if (gfld%ibmap.ne.0 .AND. gfld%ibmap.ne.254) then
258 if (.NOT. associated(gfld%bmap)) allocate(gfld%bmap(1))
259 endif
260 !
261 ! Add field to GRIB message
262 !
263 call addfield(cgrib,lcgrib,gfld%ipdtnum,gfld%ipdtmpl, &
264 gfld%ipdtlen,gfld%coord_list,gfld%num_coord, &
265 idrsnum,idrstmpl,200, &
266 gfld%fld,gfld%ngrdpts,gfld%ibmap,gfld%bmap,ierr)
267 if (ierr.ne.0) then
268 write(6,*) ' ERROR adding GRIB2 field = ',ierr
269 cycle
270 endif
271
272 enddo
273
274 if (open_grb) then ! close last GRIB2 message
275 call gribend(cgrib,lcgrib,lengrib,ierr)
276 if (ierr.ne.0) then
277 write(6,*) ' ERROR ending new GRIB2 message = ',ierr
278 if (associated(cgrib)) deallocate(cgrib)
279 call gf_free(gfld)
280 call gf_free(prevfld)
281 return
282 endif
283 open_grb=.false.
284 call wryte(ifl2,lengrib,cgrib)
285 endif
286
287 if (associated(cgrib)) deallocate(cgrib)
288 call gf_free(gfld)
289 call gf_free(prevfld)
290
291 return
292end subroutine cnv22
subroutine cnv22(ifl1, ifl2, ipack, usemiss, imiss, table_ver)
This subroutine converts every GRIB2 field in a file to another GRIB2 field, most likely one using a ...
Definition cnv22.F90:38
subroutine rdieee(rieee, a, num)
Copy array of 32-bit IEEE floating point values to local floating point representation.
Definition g2bytes.F90:637
subroutine mkieee(a, rieee, num)
Copy an array of real to an array of 32-bit IEEE floating points.
Definition g2bytes.F90:685
subroutine addgrid(cgrib, lcgrib, igds, igdstmpl, igdstmplen, ideflist, idefnum, ierr)
Add a Grid Definition Section (Section 3) to a GRIB2 message.
Definition g2create.F90:682
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.
subroutine getgb2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, unpack, k, gfld, iret)
This is a legacy version of getgb2i2().
Definition g2getgb2.F90:64
subroutine gf_free(gfld)
Free memory that was used to store array values in derived type grib_mod::gribfield.
Definition g2gf.F90:585
This Fortran module contains the declaration of derived type gribfield.
Definition gribmod.F90:10
This Fortran Module contains info on all the available GRIB Parameters, and their GRIB1 and GRIB2 cod...
Definition params.F90:12
character(len=8) function param_get_abbrev(g2disc, g2cat, g2num)
This function returns the parameter abbreviation for a given GRIB2 Discipline, Category and Parameter...
Definition params.F90:1120
Reallocate memory, preserving contents.
Definition realloc.F90:12