NCEPLIBS-g2  3.4.5
addfield.f
Go to the documentation of this file.
1 
6 
88  subroutine addfield(cgrib,lcgrib,ipdsnum,ipdstmpl,ipdstmplen,
89  & coordlist,numcoord,idrsnum,idrstmpl,
90  & idrstmplen,fld,ngrdpts,ibmap,bmap,ierr)
92  use drstemplates
93  logical :: match
94  character(len=1),intent(inout) :: cgrib(lcgrib)
95  integer,intent(in) :: ipdsnum,ipdstmpl(*)
96  integer,intent(in) :: idrsnum,numcoord,ipdstmplen,idrstmplen
97  integer,intent(in) :: lcgrib,ngrdpts,ibmap
98  real,intent(in) :: coordlist(numcoord)
99  real,target,intent(in) :: fld(ngrdpts)
100  integer,intent(out) :: ierr
101  integer,intent(inout) :: idrstmpl(*)
102  logical*1,intent(in) :: bmap(ngrdpts)
103 
104  character(len=4),parameter :: grib='GRIB',c7777='7777'
105  character(len=4):: ctemp
106  character(len=1),allocatable :: cpack(:)
107  real,pointer,dimension(:) :: pfld
108  real(4) :: coordieee(numcoord),re00
109  integer(4) :: ire00,allones
110  integer :: mappds(ipdstmplen),intbmap(ngrdpts),mapdrs(idrstmplen)
111  integer,parameter :: zero=0,one=1,four=4,five=5,six=6,seven=7
112  integer,parameter :: minsize=50000
113  integer iofst,ibeg,lencurr,len,mappdslen,mapdrslen,lpos3
114  integer width,height,ndpts
115  integer lensec3,lensec4,lensec5,lensec6,lensec7
116  logical issec3,needext,isprevbmap
117 
118  allones=z'FFFFFFFF'
119  ierr=0
120 !
121 ! Check to see if beginning of GRIB message exists
122 !
123  match=.true.
124  do i=1,4
125  if(cgrib(i) /= grib(i:i)) then
126  match=.false.
127  endif
128  enddo
129  if ( .not. match ) then
130  print *,'addfield: GRIB not found in given message.'
131  print *,'addfield: Call to routine gribcreate required',
132  & ' to initialize GRIB messge.'
133  ierr=1
134  return
135  endif
136 !
137 ! Get current length of GRIB message
138 !
139  call g2_gbytec(cgrib,lencurr,96,32)
140 !
141 ! Check to see if GRIB message is already complete
142 !
143  ctemp=cgrib(lencurr-3)//cgrib(lencurr-2)//cgrib(lencurr-1)
144  & //cgrib(lencurr)
145  if ( ctemp.eq.c7777 ) then
146  print *,'addfield: GRIB message already complete. Cannot',
147  & ' add new section.'
148  ierr=2
149  return
150  endif
151 !
152 ! Loop through all current sections of the GRIB message to
153 ! find the last section number.
154 !
155  issec3=.false.
156  isprevbmap=.false.
157  len=16 ! length of Section 0
158  do
159  ! Get number and length of next section
160  iofst=len*8
161  call g2_gbytec(cgrib,ilen,iofst,32)
162  iofst=iofst+32
163  call g2_gbytec(cgrib,isecnum,iofst,8)
164  iofst=iofst+8
165  ! Check if previous Section 3 exists and save location of
166  ! the section 3 in case needed later.
167  if (isecnum.eq.3) then
168  issec3=.true.
169  lpos3=len+1
170  lensec3=ilen
171  endif
172  ! Check if a previous defined bitmap exists
173  if (isecnum.eq.6) then
174  call g2_gbytec(cgrib,ibmprev,iofst,8)
175  iofst=iofst+8
176  if ((ibmprev.ge.0).and.(ibmprev.le.253)) isprevbmap=.true.
177  endif
178  len=len+ilen
179  ! Exit loop if last section reached
180  if ( len.eq.lencurr ) exit
181  ! If byte count for each section does not match current
182  ! total length, then there is a problem.
183  if ( len.gt.lencurr ) then
184  print *,'addfield: Section byte counts don''t add to total.'
185  print *,'addfield: Sum of section byte counts = ',len
186  print *,'addfield: Total byte count in Section 0 = ',lencurr
187  ierr=3
188  return
189  endif
190  enddo
191 !
192 ! Sections 4 through 7 can only be added after section 3 or 7.
193 !
194  if ( (isecnum.ne.3) .and. (isecnum.ne.7) ) then
195  print *,'addfield: Sections 4-7 can only be added after',
196  & ' Section 3 or 7.'
197  print *,'addfield: Section ',isecnum,' was the last found in',
198  & ' given GRIB message.'
199  ierr=4
200  return
201 !
202 ! Sections 4 through 7 can only be added if section 3 was previously defined.
203 !
204  elseif (.not.issec3) then
205  print *,'addfield: Sections 4-7 can only be added if Section',
206  & ' 3 was previously included.'
207  print *,'addfield: Section 3 was not found in',
208  & ' given GRIB message.'
209  print *,'addfield: Call to routine addgrid required',
210  & ' to specify Grid definition.'
211  ierr=6
212  return
213  endif
214 !
215 ! Add Section 4 - Product Definition Section
216 !
217  ibeg=lencurr*8 ! Calculate offset for beginning of section 4
218  iofst=ibeg+32 ! leave space for length of section
219  call g2_sbytec(cgrib,four,iofst,8) ! Store section number ( 4 )
220  iofst=iofst+8
221  call g2_sbytec(cgrib,numcoord,iofst,16) ! Store num of coordinate values
222  iofst=iofst+16
223  call g2_sbytec(cgrib,ipdsnum,iofst,16) ! Store Prod Def Template num.
224  iofst=iofst+16
225  !
226  ! Get Product Definition Template
227  !
228  call getpdstemplate(ipdsnum,mappdslen,mappds,needext,iret)
229  if (iret.ne.0) then
230  ierr=5
231  return
232  endif
233  !
234  ! Extend the Product Definition Template, if necessary.
235  ! The number of values in a specific template may vary
236  ! depending on data specified in the "static" part of the
237  ! template.
238  !
239  if ( needext ) then
240  call extpdstemplate(ipdsnum,ipdstmpl,mappdslen,mappds)
241  endif
242  !
243  ! Pack up each input value in array ipdstmpl into the
244  ! the appropriate number of octets, which are specified in
245  ! corresponding entries in array mappds.
246  !
247  do i=1,mappdslen
248  nbits=iabs(mappds(i))*8
249  if ( (mappds(i).ge.0).or.(ipdstmpl(i).ge.0) ) then
250  call g2_sbytec(cgrib,ipdstmpl(i),iofst,nbits)
251  else
252  call g2_sbytec(cgrib,one,iofst,1)
253  call g2_sbytec(cgrib,iabs(ipdstmpl(i)),iofst+1,nbits-1)
254  endif
255  iofst=iofst+nbits
256  enddo
257  !
258  ! Add Optional list of vertical coordinate values
259  ! after the Product Definition Template, if necessary.
260  !
261  if ( numcoord .ne. 0 ) then
262  call mkieee(coordlist,coordieee,numcoord)
263  call g2_sbytesc(cgrib,coordieee,iofst,32,0,numcoord)
264  iofst=iofst+(32*numcoord)
265  endif
266  !
267  ! Calculate length of section 4 and store it in octets
268  ! 1-4 of section 4.
269  !
270  lensec4=(iofst-ibeg)/8
271  call g2_sbytec(cgrib,lensec4,ibeg,32)
272 !
273 ! Pack Data using appropriate algorithm
274 !
275  !
276  ! Get Data Representation Template
277  !
278  call getdrstemplate(idrsnum,mapdrslen,mapdrs,needext,iret)
279  if (iret.ne.0) then
280  ierr=5
281  return
282  endif
283  !
284  ! contract data field, removing data at invalid grid points,
285  ! if bit-map is provided with field.
286  !
287  if ( ibmap.eq.0 .OR. ibmap.eq.254 ) then
288  allocate(pfld(max(2,ngrdpts)))
289  ndpts=0;
290  do jj=1,ngrdpts
291  intbmap(jj)=0
292  if ( bmap(jj) ) then
293  intbmap(jj)=1
294  ndpts=ndpts+1
295  pfld(ndpts)=fld(jj);
296  endif
297  enddo
298  if(ndpts==0 .and. ngrdpts>0) then
299  pfld(1)=0
300  endif
301  else
302  ndpts=ngrdpts;
303  pfld=>fld;
304  endif
305  lcpack=0
306  nsize=ndpts*4
307  if (nsize .lt. minsize) nsize=minsize
308  allocate(cpack(nsize),stat=istat)
309  if (idrsnum.eq.0) then ! Simple Packing
310  call simpack(pfld,ndpts,idrstmpl,cpack,lcpack)
311  elseif (idrsnum.eq.2.or.idrsnum.eq.3) then ! Complex Packing
312  call cmplxpack(pfld,ndpts,idrsnum,idrstmpl,cpack,lcpack)
313  elseif (idrsnum.eq.50) then ! Sperical Harmonic Simple Packing
314  call simpack(pfld(2),ndpts-1,idrstmpl,cpack,lcpack)
315  call mkieee(real(pfld(1)),re00,1) ! ensure RE(0,0) value is IEEE format
316  !call g2_gbytec(re00,idrstmpl(5),0,32)
317  ire00=transfer(re00,ire00)
318  idrstmpl(5)=ire00
319  elseif (idrsnum.eq.51) then ! Sperical Harmonic Complex Packing
320  call getpoly(cgrib(lpos3),lensec3,jj,kk,mm)
321  if (jj.ne.0 .AND. kk.ne.0 .AND. mm.ne.0) then
322  call specpack(pfld,ndpts,jj,kk,mm,idrstmpl,cpack,lcpack)
323  else
324  print *,'addfield: Cannot pack DRT 5.51.'
325  ierr=9
326  return
327  endif
328 
329  elseif (idrsnum.eq.40 .OR. idrsnum.eq.40000) then ! JPEG2000 encoding
330  if (ibmap.eq.255) then
331  call getdim(cgrib(lpos3),lensec3,width,height,iscan)
332  if ( width.eq.0 .OR. height.eq.0 ) then
333  width=ndpts
334  height=1
335  elseif ( width.eq.allones .OR. height.eq.allones ) then
336  width=ndpts
337  height=1
338  elseif ( ibits(iscan,5,1) .eq. 1) then ! Scanning mode: bit 3
339  itemp=width
340  width=height
341  height=itemp
342  endif
343  else
344  width=ndpts
345  height=1
346  endif
347  if(width<1 .or. height<1) then
348  ! Special case: bitmask off everywhere.
349  write(0,*) 'Warning: bitmask off everywhere.'
350  write(0,*) ' Pretend one point in jpcpack to avoid crash.'
351  width=1
352  height=1
353  endif
354  lcpack=nsize
355  !print *,'w,h=',width,height
356  call jpcpack(pfld,width,height,idrstmpl,cpack,lcpack)
357 
358 
359  elseif (idrsnum.eq.41 .OR. idrsnum.eq.40010) then ! PNG encoding
360  if (ibmap.eq.255) then
361  call getdim(cgrib(lpos3),lensec3,width,height,iscan)
362  if ( width.eq.0 .OR. height.eq.0 ) then
363  width=ndpts
364  height=1
365  elseif ( width.eq.allones .OR. height.eq.allones ) then
366  width=ndpts
367  height=1
368  elseif ( ibits(iscan,5,1) .eq. 1) then ! Scanning mode: bit 3
369  itemp=width
370  width=height
371  height=itemp
372  endif
373  else
374  width=ndpts
375  height=1
376  endif
377  !print *,'png size ',width,height
378  call pngpack(pfld,width,height,idrstmpl,cpack,lcpack)
379  !print *,'png packed'
380  else
381  print *,'addfield: Data Representation Template 5.',idrsnum,
382  * ' not yet implemented.'
383  ierr=7
384  return
385  endif
386  if ( ibmap.eq.0 .OR. ibmap.eq.254 ) then
387  deallocate(pfld)
388  endif
389  if ( lcpack .lt. 0 ) then
390  if( allocated(cpack) )deallocate(cpack)
391  ierr=10
392  return
393  endif
394 
395 !
396 ! Add Section 5 - Data Representation Section
397 !
398  ibeg=iofst ! Calculate offset for beginning of section 5
399  iofst=ibeg+32 ! leave space for length of section
400  call g2_sbytec(cgrib,five,iofst,8) ! Store section number ( 5 )
401  iofst=iofst+8
402  call g2_sbytec(cgrib,ndpts,iofst,32) ! Store num of actual data points
403  iofst=iofst+32
404  call g2_sbytec(cgrib,idrsnum,iofst,16) ! Store Data Repr. Template num.
405  iofst=iofst+16
406  !
407  ! Pack up each input value in array idrstmpl into the
408  ! the appropriate number of octets, which are specified in
409  ! corresponding entries in array mapdrs.
410  !
411  do i=1,mapdrslen
412  nbits=iabs(mapdrs(i))*8
413  if ( (mapdrs(i).ge.0).or.(idrstmpl(i).ge.0) ) then
414  call g2_sbytec(cgrib,idrstmpl(i),iofst,nbits)
415  else
416  call g2_sbytec(cgrib,one,iofst,1)
417  call g2_sbytec(cgrib,iabs(idrstmpl(i)),iofst+1,nbits-1)
418  endif
419  iofst=iofst+nbits
420  enddo
421  !
422  ! Calculate length of section 5 and store it in octets
423  ! 1-4 of section 5.
424  !
425  lensec5=(iofst-ibeg)/8
426  call g2_sbytec(cgrib,lensec5,ibeg,32)
427 
428 !
429 ! Add Section 6 - Bit-Map Section
430 !
431  ibeg=iofst ! Calculate offset for beginning of section 6
432  iofst=ibeg+32 ! leave space for length of section
433  call g2_sbytec(cgrib,six,iofst,8) ! Store section number ( 6 )
434  iofst=iofst+8
435  call g2_sbytec(cgrib,ibmap,iofst,8) ! Store Bit Map indicator
436  iofst=iofst+8
437  !
438  ! Store bitmap, if supplied
439  !
440  if (ibmap.eq.0) then
441  call g2_sbytesc(cgrib,intbmap,iofst,1,0,ngrdpts) ! Store BitMap
442  iofst=iofst+ngrdpts
443  endif
444  !
445  ! If specifying a previously defined bit-map, make sure
446  ! one already exists in the current GRIB message.
447  !
448  if ((ibmap.eq.254).and.(.not.isprevbmap)) then
449  print *,'addfield: Requested previously defined bitmap, ',
450  & ' but one does not exist in the current GRIB message.'
451  ierr=8
452  return
453  endif
454  !
455  ! Calculate length of section 6 and store it in octets
456  ! 1-4 of section 6. Pad to end of octect, if necessary.
457  !
458  left=8-mod(iofst,8)
459  if (left.ne.8) then
460  call g2_sbytec(cgrib,zero,iofst,left) ! Pad with zeros to fill Octet
461  iofst=iofst+left
462  endif
463  lensec6=(iofst-ibeg)/8
464  call g2_sbytec(cgrib,lensec6,ibeg,32)
465 
466 !
467 ! Add Section 7 - Data Section
468 !
469  ibeg=iofst ! Calculate offset for beginning of section 7
470  iofst=ibeg+32 ! leave space for length of section
471  call g2_sbytec(cgrib,seven,iofst,8) ! Store section number ( 7 )
472  iofst=iofst+8
473  ! Store Packed Binary Data values, if non-constant field
474  if (lcpack.ne.0) then
475  ioctet=iofst/8
476  cgrib(ioctet+1:ioctet+lcpack)=cpack(1:lcpack)
477  iofst=iofst+(8*lcpack)
478  endif
479  !
480  ! Calculate length of section 7 and store it in octets
481  ! 1-4 of section 7.
482  !
483  lensec7=(iofst-ibeg)/8
484  call g2_sbytec(cgrib,lensec7,ibeg,32)
485 
486  if( allocated(cpack) )deallocate(cpack)
487 !
488 ! Update current byte total of message in Section 0
489 !
490  newlen=lencurr+lensec4+lensec5+lensec6+lensec7
491  call g2_sbytec(cgrib,newlen,96,32)
492 
493  return
494  end
pdstemplates::extpdstemplate
subroutine extpdstemplate(number, list, nummap, map)
This subroutine generates the remaining octet map for a given Product Definition Template,...
Definition: pdstemplates.f:459
drstemplates
This Fortran Module contains info on all the available GRIB2 Data Representation Templates used in Se...
Definition: drstemplates.f:37
drstemplates::getdrstemplate
subroutine getdrstemplate(number, nummap, map, needext, iret)
This subroutine returns DRS template information for a .
Definition: drstemplates.f:165
specpack
subroutine specpack(fld, ndpts, JJ, KK, MM, idrstmpl, cpack, lcpack)
This subroutine packs a spectral data field using the complex packing algorithm for spherical harmoni...
Definition: specpack.f:23
g2_sbytesc
subroutine g2_sbytesc(OUT, IN, ISKIP, NBYTE, NSKIP, N)
This subrountine is to put arbitrary size values into a packed bit string, taking the low order bits ...
Definition: g2_gbytesc.f:115
pdstemplates::getpdstemplate
subroutine getpdstemplate(number, nummap, map, needext, iret)
This subroutine returns PDS template information for a specified Product Definition Template 4....
Definition: pdstemplates.f:412
addfield
subroutine addfield(cgrib, lcgrib, ipdsnum, ipdstmpl, ipdstmplen, coordlist, numcoord, idrsnum, idrstmpl, idrstmplen, fld, ngrdpts, ibmap, bmap, ierr)
This subroutine packs up Sections 4 through 7 for a given field and adds them to a GRIB2 message.
Definition: addfield.f:91
pngpack
subroutine pngpack(fld, width, height, idrstmpl, cpack, lcpack)
This subroutine packs up a data field into PNG image format.
Definition: pngpack.f:32
simpack
subroutine simpack(fld, ndpts, idrstmpl, cpack, lcpack)
This subroutine packs up a data field using a simple packing algorithm as defined in the GRIB2 docume...
Definition: simpack.f:34
g2_sbytec
subroutine g2_sbytec(OUT, IN, ISKIP, NBYTE)
This subrountine is to put arbitrary size values into a packed bit string, taking the low order bits ...
Definition: g2_gbytesc.f:39
cmplxpack
subroutine cmplxpack(fld, ndpts, idrsnum, idrstmpl, cpack, lcpack)
This subroutine packs up a data field using a complex packing algorithm as defined in the GRIB2 docum...
Definition: cmplxpack.f:32
getdim
subroutine getdim(csec3, lcsec3, width, height, iscan)
This subroutine returns the dimensions and scanning mode of a grid definition packed in GRIB2 Grid De...
Definition: getdim.f:22
pdstemplates
This Fortran Module contains info on all the available GRIB2 Product Definition Templates used in Sec...
Definition: pdstemplates.f:54
mkieee
subroutine mkieee(a, rieee, num)
This subroutine stores a list of real values in 32-bit IEEE floating point format.
Definition: mkieee.f:17
jpcpack
subroutine jpcpack(fld, width, height, idrstmpl, cpack, lcpack)
This subroutine packs up a data field into a JPEG2000 code stream.
Definition: jpcpack.f:39
g2_gbytec
subroutine g2_gbytec(IN, IOUT, ISKIP, NBYTE)
This subrountine is to extract arbitrary size values from a packed bit string, right justifying each ...
Definition: g2_gbytesc.f:20
getpoly
subroutine getpoly(csec3, lcsec3, jj, kk, mm)
This subroutine returns the J, K, and M pentagonal resolution parameters specified in a GRIB Grid Def...
Definition: getpoly.f:24