NCEPLIBS-g2 4.0.0
Loading...
Searching...
No Matches
copygb2.F90
Go to the documentation of this file.
1
4
161PROGRAM copygb2
162 CHARACTER*256 carg,cg1,cx1,cgb,cxb,cgm,cxm,cg2,cnl
163 INTEGER karg(100)
164 INTEGER kgdti(200),ipopt(20),jpdt(200),jpdsb(200),iuv(100)
165 DATA igdtn/-1/,kgdti/200*0/
166 DATA ip/0/,ipopt/20*-1/
167 DATA jpdtn/-1/,jpdt/200*-9999/,jpdsb/200*-1/
168 DATA iuv/514,99*0/,nuv/1/
169 DATA lwg/0/,lapp/0/,lxx/0/,lx/1/,kz1/-1/,kz2/-2/
170 DATA jb/0/,jbk/0/,lab/1/,ab/-1.e30/,lam/0/,am/0./
171 DATA cgb/' '/,cxb/' '/,cgm/' '/,cxm/' '/,cnl/' '/
172 INTEGER ids(255),ibs(255),nbs(255)
173 NAMELIST/nlcopygb/ ids,ibs,nbs
174 DATA ids/255*-9999/,ibs/255*-9999/,nbs/255*-9999/
175
176 ! PARSE COMMAND LINE OPTIONS
177 narg=iargc()
178 iarg=1
179 lstopt=0
180 DO WHILE(iarg.LE.narg.AND.lstopt.EQ.0)
181 CALL getarg(iarg,carg)
182 larg=len_trim(carg)
183 iarg=iarg+1
184 IF(carg(1:1).NE.'-') THEN
185 lstopt=1
186 iarg=iarg-1
187 ELSEIF(larg.EQ.1) THEN
188 CALL errmsg('copygb2: invalid option -')
189 CALL eusage
190 CALL errexit(1)
191 ELSE
192 l=2
193 DO WHILE(l.LE.larg)
194 IF(carg(l:l).EQ.'-') THEN
195 lstopt=1
196 ELSEIF(carg(l:l).EQ.'a') THEN
197 lapp=1
198 ELSEIF(carg(l:l).EQ.'A') THEN
199 IF(l.EQ.larg) THEN
200 l=0
201 CALL getarg(iarg,carg)
202 larg=len_trim(carg)
203 iarg=iarg+1
204 ENDIF
205 IF(carg(l+1:l+1).EQ.'>') THEN
206 lab=1
207 l=l+1
208 ELSEIF(carg(l+1:l+1).EQ.'<') THEN
209 lab=-1
210 l=l+1
211 ELSE
212 CALL errmsg('copygb2: invalid threshold '// &
213 carg(l+1:larg))
214 CALL eusage
215 CALL errexit(1)
216 ENDIF
217 CALL fparser(carg(l+1:larg),1,ab)
218 l=larg
219 call errmsg('Option -A Ignored...Not yet implemented.')
220 lab=1 ! default value, since -A option not yet implemented.
221 ab=-1.e30 ! default value, since -A option not yet implemented.
222 ELSEIF(carg(l:l).EQ.'B') THEN
223 IF(l.EQ.larg) THEN
224 l=0
225 CALL getarg(iarg,carg)
226 larg=len_trim(carg)
227 iarg=iarg+1
228 ENDIF
229 lcgb=larg-l
230 cgb=carg(l+1:larg)
231 l=larg
232 call errmsg('Option -B Ignored...Not yet implemented.')
233 lcgb=1 ! default value, since -B option not yet implemented.
234 cgb=' ' ! default value, since -B option not yet implemented.
235 ELSEIF(carg(l:l).EQ.'b') THEN
236 IF(l.EQ.larg) THEN
237 l=0
238 CALL getarg(iarg,carg)
239 larg=len_trim(carg)
240 iarg=iarg+1
241 ENDIF
242 lcxb=larg-l
243 cxb=carg(l+1:larg)
244 l=larg
245 call errmsg('Option -b Ignored...Not yet implemented.')
246 lcxb=1 ! default value, since -B option not yet implemented.
247 cxb=' ' ! default value, since -B option not yet implemented.
248 ELSEIF(carg(l:l).EQ.'g') THEN
249 IF(l.EQ.larg) THEN
250 l=0
251 CALL getarg(iarg,carg)
252 larg=len_trim(carg)
253 iarg=iarg+1
254 ENDIF
255 karg(1)=igdtn
256 karg(2:100)=0
257 CALL fparsei(carg(l+1:larg),100,karg)
258 igdtn=karg(1)
259 IF(igdtn.GE.0.AND.igdtn.LE.65534) THEN
260 kgdti(1:99)=karg(2:100)
261 ENDIF
262 IF(igdtn.LT.-5.OR.igdtn.EQ.-2.OR. &
263 igdtn.EQ.-3.OR.igdtn.GT.65534) THEN
264 CALL errmsg('copygb2: invalid output grid '// &
265 carg(l+1:larg))
266 CALL eusage
267 CALL errexit(1)
268 ENDIF
269 IF ( igdtn.GE.0 ) THEN
270 mi=numpts(igdtn,kgdti)
271 IF(mi.LE.0) THEN
272 CALL errmsg('copygb2: unsupported output grid '// &
273 carg(l+1:larg))
274 CALL eusage
275 CALL errexit(1)
276 ENDIF
277 ENDIF
278 l=larg
279 ELSEIF(carg(l:l).EQ.'i') THEN
280 IF(l.EQ.larg) THEN
281 l=0
282 CALL getarg(iarg,carg)
283 larg=len_trim(carg)
284 iarg=iarg+1
285 ENDIF
286 karg(1)=ip
287 karg(2:21)=ipopt
288 CALL fparsei(carg(l+1:larg),21,karg)
289 ip=karg(1)
290 ipopt=karg(2:21)
291 l=larg
292 ELSEIF(carg(l:l).EQ.'K') THEN
293 IF(l.EQ.larg) THEN
294 l=0
295 CALL getarg(iarg,carg)
296 larg=len_trim(carg)
297 iarg=iarg+1
298 ENDIF
299 jbk=1
300 CALL fparsei(carg(l+1:larg),100,jpdsb)
301 IF(jpdsb(5).EQ.0) THEN
302 CALL errmsg('copygb2: invalid PDS parms '// &
303 carg(l+1:larg))
304 CALL eusage
305 CALL errexit(1)
306 ENDIF
307 l=larg
308 call errmsg('Option -K Ignored...Not yet implemented.')
309 jbk=0 ! default value, since -K option not yet implemented.
310 jpdsb=-1 ! default value, since -K option not yet implemented.
311 ELSEIF(carg(l:l).EQ.'k') THEN
312 IF(l.EQ.larg) THEN
313 l=0
314 CALL getarg(iarg,carg)
315 larg=len_trim(carg)
316 iarg=iarg+1
317 ENDIF
318 karg(1)=jpdtn
319 karg(2:100)=jpdt(1:99)
320 CALL fparsei(carg(l+1:larg),100,karg)
321 jpdtn=karg(1)
322 jpdt(1:99)=karg(2:100)
323 IF(jpdtn.LT.-1 .OR. jpdtn.GE.65535) THEN
324 CALL errmsg('copygb2: invalid PDT parms '// &
325 carg(l+1:larg))
326 CALL eusage
327 CALL errexit(1)
328 ENDIF
329 l=larg
330 ELSEIF(carg(l:l).EQ.'M') THEN
331 IF(l.EQ.larg) THEN
332 l=0
333 CALL getarg(iarg,carg)
334 larg=len_trim(carg)
335 iarg=iarg+1
336 ENDIF
337 IF(carg(l+1:l+1).EQ.'#') THEN
338 l=l+1
339 CALL fparser(carg(l+1:larg),1,am)
340 lam=1
341 ELSE
342 lcgm=larg-l
343 cgm=carg(l+1:larg)
344 lam=5
345 ENDIF
346 l=larg
347 ELSEIF(carg(l:l).EQ.'m') THEN
348 IF(l.EQ.larg) THEN
349 l=0
350 CALL getarg(iarg,carg)
351 larg=len_trim(carg)
352 iarg=iarg+1
353 ENDIF
354 lcxm=larg-l
355 cxm=carg(l+1:larg)
356 l=larg
357 ELSEIF(carg(l:l).EQ.'N') THEN
358 IF(l.EQ.larg) THEN
359 l=0
360 CALL getarg(iarg,carg)
361 larg=len_trim(carg)
362 iarg=iarg+1
363 ENDIF
364 lcnl=larg-l
365 cnl=carg(l+1:larg)
366 l=larg
367 call errmsg('Option -N Ignored...Not yet implemented.')
368 lcnl=1 ! default value, since -N option not yet implemented.
369 cnl=' ' ! default value, since -N option not yet implemented.
370 ELSEIF(carg(l:l).EQ.'v') THEN
371 IF(l.EQ.larg) THEN
372 l=0
373 CALL getarg(iarg,carg)
374 larg=len_trim(carg)
375 iarg=iarg+1
376 ENDIF
377 CALL fparsei(carg(l+1:larg),100,iuv)
378 nuv=1
379 DO juv=2,100
380 IF(iuv(juv).NE.0) nuv=juv
381 ENDDO
382 l=larg
383 ELSEIF(carg(l:l).EQ.'x') THEN
384 lx=0
385 ELSEIF(carg(l:l).EQ.'X') THEN
386 lxx=1
387 ELSE
388 CALL errmsg('copygb2: invalid option '//carg(l:l))
389 CALL eusage
390 CALL errexit(1)
391 ENDIF
392 l=l+1
393 ENDDO
394 ENDIF
395 ENDDO
396
397 ! PARSE COMMAND LINE POSITIONAL ARGUMENTS
398 nxarg=lx+2
399 IF(narg-iarg+1.NE.nxarg) THEN
400 CALL errmsg('copygb2: incorrect number of arguments')
401 CALL eusage
402 CALL errexit(nxarg)
403 ENDIF
404 CALL getarg(iarg,cg1)
405 lcg1=len_trim(cg1)
406 iarg=iarg+1
407 lg1=11
408 CALL baopenr(lg1,cg1(1:lcg1),iretba)
409 IF(iretba.NE.0) THEN
410 CALL errmsg('copygb2: error accessing file '//cg1(1:lcg1))
411 CALL errexit(8)
412 ENDIF
413 IF(lx.GT.0) THEN
414 CALL getarg(iarg,cx1)
415 lcx1=len_trim(cx1)
416 iarg=iarg+1
417 lx1=31
418 CALL baopenr(lx1,cx1(1:lcx1),iretba)
419 IF(iretba.NE.0) THEN
420 CALL errmsg('copygb2: error accessing file '//cx1(1:lcx1))
421 CALL errexit(8)
422 ENDIF
423 ELSE
424 lx1=0
425 ENDIF
426 CALL getarg(iarg,cg2)
427 lcg2=len_trim(cg2)
428 iarg=iarg+1
429 IF(cg2(1:lcg2).EQ.'-') THEN
430 IF(lxx.GT.0) THEN
431 CALL errmsg('copygb2: piping incompatible with the X option')
432 CALL errexit(1)
433 ENDIF
434 lg2=6
435 ELSE
436 lg2=51
437 IF(lapp.EQ.0) THEN
438 CALL baopenwt(lg2,cg2(1:lcg2),iretba)
439 ELSE
440 CALL baopenwa(lg2,cg2(1:lcg2),iretba)
441 ENDIF
442 IF(iretba.NE.0) THEN
443 CALL errmsg('copygb2: error accessing file '//cg2(1:lcg2))
444 CALL errexit(8)
445 ENDIF
446 ENDIF
447
448 ! OPEN MAP FILE
449 IF(cgb.NE.' ') THEN
450 IF(cgb(1:2).EQ.'-1') THEN
451 IF(jpdsb(5).EQ.-1) THEN
452 jb=1
453 ELSE
454 jb=4
455 lgb=lg1
456 lxb=lx1
457 ENDIF
458 ELSE
459 jb=4
460 lgb=14
461 CALL baopenr(lgb,cgb(1:lcgb),iretba)
462 IF(iretba.NE.0) THEN
463 CALL errmsg('copygb2: error accessing file '//cgb(1:lcgb))
464 CALL errexit(8)
465 ENDIF
466 IF(cxb(1:1).NE.' ') THEN
467 lxb=34
468 CALL baopenr(lxb,cxb(1:lcxb),iretba)
469 IF(iretba.NE.0) THEN
470 CALL errmsg('copygb2: error accessing file '//cxb(1:lcxb))
471 CALL errexit(8)
472 ENDIF
473 ELSE
474 lxb=0
475 ENDIF
476 ENDIF
477 ENDIF
478
479 ! OPEN MERGE FILE
480 IF(cgm.NE.' ') THEN
481 lam=5
482 lgm=15
483 CALL baopenr(lgm,cgm(1:lcgm),iretba)
484 IF(iretba.NE.0) THEN
485 CALL errmsg('copygb2: error accessing file '//cgm(1:lcgm))
486 CALL errexit(8)
487 ENDIF
488 IF(cxm(1:1).NE.' ') THEN
489 lxm=35
490 CALL baopenr(lxm,cxm(1:lcxm),iretba)
491 IF(iretba.NE.0) THEN
492 CALL errmsg('copygb2: error accessing file '//cxm(1:lcxm))
493 CALL errexit(8)
494 ENDIF
495 ELSE
496 lxm=0
497 ENDIF
498 ENDIF
499
500 ! OPEN AND READ NAMELIST FILE
501 IF(cnl.NE.' ') THEN
502 lnl=2
503 OPEN(lnl,file=cnl(1:lcnl),status='OLD',iostat=iret)
504 IF(iret.NE.0) THEN
505 CALL errmsg('copygb2: error accessing file '//cnl(1:lcnl))
506 CALL errexit(8)
507 ENDIF
508 READ(lnl,nlcopygb,iostat=iret)
509 IF(iret.NE.0) THEN
510 CALL errmsg('copygb2: error reading namelist from file '// &
511 cnl(1:lcnl))
512 CALL errexit(8)
513 ENDIF
514 ENDIF
515
516 ! GO
517 IF(lxx.GT.0) THEN
518 CALL w3tagb('COPYGB2 ',1998,0295,0047,'NP23 ')
519 ENDIF
520 CALL cpgb(lg1,lx1,lgb,lxb,lgm,lxm,lg2, &
521 igdtn,kgdti,ip,ipopt,jpdtn,jpdt,nuv,iuv, &
522 jpdsb,jb,jbk,lab,ab,lam,am,lxx,lwg, &
523 ids,ibs,nbs)
524 IF(lxx.GT.0) THEN
525 CALL w3tage('COPYGB2 ')
526 ENDIF
527
528CONTAINS
529
533 SUBROUTINE eusage
534 CALL errmsg('Usage: copygb2'// &
535 ' [-g "kgdtn [kgdt]"] [-i "ip [ipopts]"]')
536 CALL errmsg(' '// &
537 ' [-k "kpdtn kpdt"] [-v "uparms"]')
538 CALL errmsg(' '// &
539 ' [-B mapgrib [-b mapindex] [-A "<> mapthreshold"]'// &
540 ' [-K "mapkpds"]]')
541 CALL errmsg(' '// &
542 ' [-M "mask"/mergegrib [-m mergeindex]] [-X] [-a]'// &
543 ' [-N namelist]')
544 CALL errmsg(' then either:')
545 CALL errmsg(' '// &
546 ' grib2in index1 grib2out')
547 CALL errmsg(' or:')
548 CALL errmsg(' '// &
549 ' -x grib2in grib2out')
550
551 END SUBROUTINE eusage
552
584 SUBROUTINE cpgb(LG1,LX1,LGB,LXB,LGM,LXM,LG2, &
585 IGDTN,KGDTI,IP,IPOPT,JPDTN,JPDT,NUV,IUV, &
586 JPDSB,JB,JBK,LAB,AB,LAM,AM,LXX,LWG, &
587 IDS,IBS,NBS)
588 USE grib_mod
589
590 parameter(mbuf=256*1024)
591 CHARACTER cbufb(mbuf)
592 INTEGER jids(200),jpdt(200),jgdt(200)
593 INTEGER jjpdt(200)
594 INTEGER jpdsb(100),iuv(100)
595 INTEGER kgdti(200)
596 INTEGER ipopt(20)
597 INTEGER ids(255),ibs(255),nbs(255)
598 !INTEGER JPDS(200)
599 INTEGER jgds(200),jens(5)
600 INTEGER kpdsb(200),kgdsb(200),kensb(5)
601 !INTEGER KPDSM(200),KGDSM(200),KENSM(5)
602 !CHARACTER*80 CIN
603 LOGICAL unpack
604 TYPE( gribfield ) :: gfld1,gfldm
605
606 mm = 0
607 ! READ GRIB HEADERS
608 IF(lxx.GT.0) CALL instrument(6,kall0,ttot0,tmin0,tmax0)
609 IF(jb.EQ.4) THEN
610 jgds=-1
611 jens=-1
612 krb=-1
613 kpdsb=0
614 kgdsb=0
615 CALL getgbemh(lgb,lxb,krb,jpdsb,jgds,jens, &
616 mbuf,cbufb,nlenb,nnumb,mnumb, &
617 kb,mb,krbx,kpdsb,kgdsb,kensb,iret)
618 IF(iret.EQ.0.AND.mb.LE.0) iret=255
619 IF(lxx.GT.0) THEN
620 IF(iret.EQ.99) THEN
621 print *,'copygb2 map field not found'
622 ELSEIF(iret.NE.0) THEN
623 print *,'copygb2 map field retrieval error code ',iret
624 ENDIF
625 ENDIF
626 ELSE
627 mb=1
628 iret=0
629 ENDIF
630 IF(iret.EQ.0) THEN
631 kr1=0
632 !IF(LWG.EQ.1) THEN
633 ! READ (*,*,IOSTAT=IRET) CIN
634 ! IF(IRET.EQ.0) THEN
635 ! NDEL=SCAN(CIN,":")
636 ! IF(NDEL.GT.0) CIN=CIN(:NDEL-1)
637 ! READ(CIN,*) KR1
638 ! KR1=-KR1
639 ! ENDIF
640 !ENDIF
641 IF(iret.EQ.0) THEN
642 jdisc=-1
643 jids=-9999
644 jgdtn=-1
645 jgdt=-9999
646 unpack=.false.
647 CALL getgb2(lg1,lx1,kr1,jdisc,jids,jpdtn,jpdt,jgdtn, &
648 jgdt,unpack,kr1x,gfld1,iret)
649 ! CALL GETGBEMH(LG1,LX1,KR1,JPDS1,JGDS,JENS, &
650 ! MBUF,CBUF1,NLEN1,NNUM1,MNUM1, &
651 ! K1,M1,KR1X,KPDS1,KGDS1,KENS1,IRET)
652 IF(iret.EQ.0.AND.gfld1%NDPTS.LE.0) iret=255
653 m1=gfld1%NGRDPTS
654 kr1=kr1x
655 IF(lxx.GT.0) THEN
656 IF(iret.EQ.99) THEN
657 print *,'copygb2 field not found'
658 ELSEIF(iret.NE.0) THEN
659 print *,'copygb2 header retrieval error code ',iret
660 ENDIF
661 ENDIF
662 ENDIF
663 ENDIF
664
665 ! LOOP UNTIL DONE
666 no=0
667 DO WHILE(iret.EQ.0)
668 IF(lam.EQ.5) THEN
669 jdisc=gfld1%DISCIPLINE
670 jids=-9999
671 jjpdtn=gfld1%IPDTNUM
672 jjpdt=-9999
673 jjpdt(1:2)=gfld1%IPDTMPL(1:2)
674 jjpdt(10:15)=gfld1%IPDTMPL(10:15)
675 jgdtn=-1
676 jgdt=-9999
677 unpack=.false.
678 krm=0
679 CALL getgb2(lgm,lxm,krm,jdisc,jids,jjpdtn,jjpdt, &
680 jgdtn,jgdt,unpack,krmx,gfldm,iret)
681 ! CALL GETGBEMH(LGM,LXM,KRM,JPDS,JGDS,JENS, &
682 ! MBUF,CBUFM,NLENM,NNUMM,MNUMM, &
683 ! KM,MM,KRMX,KPDSM,KGDSM,KENSM,IRET)
684 mm=gfldm%NGRDPTS
685 IF(iret.EQ.0.AND.mm.LE.0) iret=255
686 IF(iret.NE.0) THEN
687 mm=0
688 gfldm%IGDTNUM=-1
689 iret=0
690 ENDIF
691 ENDIF
692 IF(igdtn.EQ.-1) THEN
693 igdtn=gfld1%IGDTNUM
694 kgdti(1:gfld1%IGDTLEN) = gfld1%IGDTMPL(1:gfld1%IGDTLEN)
695 mi=m1
696 ELSEIF(igdtn.EQ.-4.AND.jb.EQ.4) THEN
697 igi=kpdsb(3)
698 !IGDTN=
699 !KGDSI=KGDSB
700 !KGDTI=
701 mi=mb
702 ELSEIF(igdtn.EQ.-5.AND.lam.EQ.5) THEN
703 igdtn=gfldm%IGDTNUM
704 kgdti(1:gfldm%IGDTLEN) = gfldm%IGDTMPL(1:gfldm%IGDTLEN)
705 mi=mm
706 ELSE
707 mi=numpts(igdtn,kgdti)
708 ENDIF
709 IF(lxx.GT.0) CALL instrument(1,kall1,ttot1,tmin1,tmax1)
710 IF(igdtn.GE.0.AND.igdtn.LE.65534) THEN
711 mf=max(m1,mb,mm)
712 CALL cpgb1(lg1,lx1,m1, &
713 mbuf,mf,mi, &
714 igdtn,kgdti,ip,ipopt,jpdtn,jpdt,nuv,iuv, &
715 jpdsb,jb,jbk,lab,ab,lam,am, &
716 ids,ibs,nbs, &
717 lgb,lxb,mb,cbufb,nlenb,nnumb,mnumb, &
718 lgm,lxm,mm, &
719 lg2,lxx,kr1-1,no,iret1)
720 ENDIF
721 IF(lam.EQ.5) THEN ! clean-up
722 CALL gf_free(gfldm)
723 ENDIF
724 !IF(LWG.EQ.1) THEN
725 ! READ (*,*,IOSTAT=IRET) CIN
726 ! IF(IRET.EQ.0) THEN
727 ! NDEL=SCAN(CIN,":")
728 ! IF(NDEL.GT.0) CIN=CIN(:NDEL-1)
729 ! READ(CIN,*) KR1
730 ! KR1=KR1-1
731 ! ENDIF
732 !ENDIF
733 IF(iret.EQ.0) THEN
734 CALL gf_free(gfld1)
735 jdisc=-1
736 jids=-9999
737 jgdtn=-1
738 jgdt=-9999
739 unpack=.false.
740 CALL getgb2(lg1,lx1,kr1,jdisc,jids,jpdtn,jpdt,jgdtn, &
741 jgdt,unpack,kr1x,gfld1,iret)
742 ! CALL GETGBEMH(LG1,LX1,KR1,JPDS1,JGDS,JENS, &
743 ! MBUF,CBUF1,NLEN1,NNUM1,MNUM1, &
744 ! K1,M1,KR1X,KPDS1,KGDS1,KENS1,IRET)
745 IF(iret.EQ.0.AND.gfld1%NDPTS.LE.0) iret=255
746 m1=gfld1%NGRDPTS
747 kr1=kr1x
748 IF(lxx.GT.0) THEN
749 IF(iret.NE.0.AND.iret.NE.99) THEN
750 print *,'copygb2 header retrieval error code ',iret
751 ENDIF
752 ENDIF
753 ENDIF
754 ENDDO
755
756 IF(lxx.GT.0) THEN
757 print *,'copygb2 wrote ',no,' total records'
758 CALL instrument(1,kall1,ttot1,tmin1,tmax1)
759 print *,'Instrumentation Report'
760 print '(F10.3," seconds spent searching headers")',ttot1
761 CALL instrument(-2,kall2,ttot2,tmin2,tmax2)
762 print '(F10.3," seconds spent reading and unpacking")',ttot2
763 CALL instrument(-3,kall3,ttot3,tmin3,tmax3)
764 print '(F10.3," seconds spent manipulating masks")',ttot3
765 CALL instrument(-4,kall4,ttot4,tmin4,tmax4)
766 print '(F10.3," seconds spent interpolating or copying")',ttot4
767 CALL instrument(-5,kall5,ttot5,tmin5,tmax5)
768 print '(F10.3," seconds spent merging")',ttot5
769 CALL instrument(-6,kall6,ttot6,tmin6,tmax6)
770 print '(F10.3," seconds spent packing and writing")',ttot6
771 ttott=ttot1+ttot2+ttot3+ttot4+ttot5+ttot6
772 print '(F10.3," total seconds spent in copygb2")',ttott
773 ENDIF
774
775 END SUBROUTINE cpgb
776
820 SUBROUTINE cpgb1(LG1,LX1,M1, &
821 MBUF,MF,MI, &
822 IGDTN,KGDTI,IP,IPOPT,JPDTN,JPDT,NUV,IUV, &
823 JPDSB,JB,JBK,LAB,AB,LAM,AM, &
824 IDS,IBS,NBS, &
825 LGB,LXB,MB,CBUFB,NLENB,NNUMB,MNUMB, &
826 LGM,LXM,MM, &
827 LG2,LXX,KS1,NO,IRET)
828 USE grib_mod
829 USE gridtemplates
830
831 CHARACTER cbufb(mbuf)
832 INTEGER jpdsb(100),iuv(100)
833 INTEGER jids(200),jpdt(100),jgdt(200)
834 INTEGER jjpdt(200)
835 INTEGER,TARGET :: kgdti(200)
836 INTEGER ipopt(20)
837 INTEGER ids(255),ibs(255),nbs(255)
838 !INTEGER JPDS(200)
839 INTEGER jgds(200),jens(5)
840 !INTEGER KPDS1(200),KGDS1(200)
841 !INTEGER KENS1(5)
842 INTEGER kpdsb(200),kgdsb(200),kensb(5)
843 !INTEGER KPDSM(200),KGDSM(200),KENSM(5)
844 INTEGER,POINTER :: tmpptr(:)
845 LOGICAL*1 lr(mf)
846 LOGICAL*1,POINTER :: l1i(:),lbi(:)
847 LOGICAL unpack
848 REAL fr(mf)
849 REAL gr(mf)
850 REAL,POINTER :: f1i(:)
851 REAL,ALLOCATABLE,TARGET :: g1i(:)
852 REAL,POINTER :: fbi(:),gbi(:)
853 TYPE( gribfield ) :: gfld1,gfldv,gfldm,gfldmv
854 INTEGER isdummy,iadummy(200)
855
856 ! GET FIELD FROM FILE 1
857 jdisc=-1
858 jids=-9999
859 jgdtn=-1
860 jgdt=-9999
861 unpack=.true.
862 CALL getgb2(lg1,lx1,ks1,jdisc,jids,jpdtn,jpdt,jgdtn, &
863 jgdt,unpack,kr1,gfld1,iret)
864 k1=gfld1%NGRDPTS
865 ! CALL GETGBEM(LG1,LX1,M1,KS1,JPDS1,JGDS,JENS, &
866 ! MBUF,CBUF1,NLEN1,NNUM1,MNUM1, &
867 ! K1,KR1,KPDS1,KGDS1,KENS1,LR,FR,IRET)
868 iv=0
869 krv=0
870 IF(iret.EQ.0) THEN
871 juv=1
872 nparm=(65536*gfld1%DISCIPLINE) + (256*gfld1%IPDTMPL(1)) + &
873 gfld1%IPDTMPL(2)
874 DO WHILE(juv.LE.nuv.AND.nparm.NE.iuv(juv).AND. &
875 nparm.NE.iuv(juv)+1)
876 juv=juv+1
877 ENDDO
878 IF(juv.LE.nuv.AND.nparm.EQ.iuv(juv)) THEN
879 iv=1
880 gfld1%IPDTMPL(2)=gfld1%IPDTMPL(2)+1
881 CALL getgb2(lg1,lx1,krv,gfld1%DISCIPLINE,gfld1%IDSECT, &
882 gfld1%IPDTNUM,gfld1%IPDTMPL,gfld1%IGDTNUM, &
883 gfld1%IGDTMPL,unpack,krvx,gfldv,iret)
884 ! CALL GETGBEM(LG1,LX1,M1,KRV,JPDS,JGDS,JENS, &
885 ! MBUF,CBUF1,NLEN1,NNUM1,MNUM1, &
886 ! K1,KRVX,KPDS1,KGDS1,KENS1,LR,GR,IRET)
887 krv=krvx
888 gfld1%IPDTMPL(2)=gfld1%IPDTMPL(2)-1
889 ELSEIF(juv.LE.nuv.AND.nparm.EQ.iuv(juv)+1) THEN
890 iret=-1
891 ENDIF
892 ENDIF
893 IF(lxx.GT.0) THEN
894 IF(iret.EQ.-1) THEN
895 print *,'copygb2 skipping 2nd vector component field'
896 ELSEIF(iret.NE.0) THEN
897 print *,'copygb2 data retrieval error code ',iret
898 ELSEIF(krv.EQ.0) THEN
899 print *,'copygb2 read scalar field from record ',kr1
900 print *,' ...PDT 4.',gfld1%IPDTNUM,'=', &
901 (gfld1%IPDTMPL(i),i=1,gfld1%IPDTLEN)
902 ELSE
903 print *,'copygb2 read vector field from records ',kr1,krv
904 print *,' ...PDT 4.',gfld1%IPDTNUM,'=', &
905 (gfld1%IPDTMPL(i),i=1,gfld1%IPDTLEN)
906 print *,' ...PDT 4.',gfldv%IPDTNUM,'=', &
907 (gfldv%IPDTMPL(i),i=1,gfldv%IPDTLEN)
908 ENDIF
909 CALL instrument(2,kall2,ttot2,tmin2,tmax2)
910 ENDIF
911
912 ! INVOKE MAP MASK BEFORE INTERPOLATION
913 IF(iret.EQ.0.AND.jbk.EQ.1.AND.jb.EQ.1) THEN
914 DO i=1,k1
915 IF(lr(i)) THEN
916 IF((lab.EQ.1.AND.fr(i).LE.ab).OR. &
917 (lab.EQ.-1.AND.fr(i).GE.ab)) THEN
918 ib1=1
919 lr(i)=.false.
920 ENDIF
921 ENDIF
922 ENDDO
923 IF(lxx.GT.0) THEN
924 print *,' applied pre-interpolation map mask'
925 CALL instrument(3,kall3,ttot3,tmin3,tmax3)
926 ENDIF
927 ENDIF
928
929 ! INTERPOLATE FIELD
930 IF(iret.EQ.0) THEN
931 ALLOCATE(l1i(mi))
932 ALLOCATE(f1i(mi))
933 ALLOCATE(g1i(mi))
934 IF ( gfld1%IBMAP.EQ.0 .OR. gfld1%IBMAP.EQ.254 ) THEN
935 ib1=1
936 ELSE
937 ib1=0
938 ALLOCATE( gfld1%BMAP(k1) ) ! dummy array
939 ENDIF
940 IF ( .NOT. ASSOCIATED(gfld1%LIST_OPT)) &
941 ALLOCATE(gfld1%LIST_OPT(1))
942 IF ( .NOT. ASSOCIATED(gfldv%FLD) ) ALLOCATE(gfldv%FLD(k1))
943 CALL intgrib2(iv,ip,ipopt,gfld1%IGDTNUM,gfld1%IGDTMPL, &
944 gfld1%NUM_OPT,gfld1%LIST_OPT, &
945 k1,ib1,gfld1%BMAP, &
946 gfld1%FLD,gfldv%FLD,igdtn,kgdti,mi, &
947 ib1i,l1i,f1i,g1i,iret)
948 IF(lxx.GT.0) THEN
949 IF(iret.EQ.0) THEN
950 print *,' interpolated to grid GDT 3.',igdtn,'=', &
951 (kgdti(j),j=1,getgdtlen(igdtn))
952 ELSEIF(iret.GT.0) THEN
953 print *,' interpolation error code ',iret
954 ENDIF
955 CALL instrument(4,kall4,ttot4,tmin4,tmax4)
956 ENDIF
957 IF(iret.EQ.-1) iret=0
958 ENDIF
959
960 ! GET MAP FIELD
961 IF(iret.EQ.0.AND.jb.EQ.4) THEN
962 krb=0
963 jgds=-1
964 jens=-1
965 CALL getgbem(lgb,lxb,mb,krb,jpdsb,jgds,jens, &
966 mbuf,cbufb,nlenb,nnumb,mnumb, &
967 kb,krbx,kpdsb,kgdsb,kensb,lr,fr,iret)
968 IF(lxx.GT.0) THEN
969 IF(iret.EQ.0) THEN
970 print *,' map field retrieved'
971 print *,' ...KPDS(1:24)=',(kpdsb(i),i=1,24)
972 ELSEIF(iret.EQ.99) THEN
973 print *,' map field not found'
974 ELSE
975 print *,' map field retrieval error code ',iret
976 ENDIF
977 ENDIF
978
979 ! INTERPOLATE MAP FIELD
980 IF(iret.EQ.0) THEN
981 ibb=mod(kpdsb(4)/64,2)
982 CALL intgrib2(0,ip,ipopt,isdummy,iadummy,SIZE(kgdsb),kgdsb, &
983 kb,ibb,lr,fr,gr,isdummy,kgdti,mi, &
984 ibbi,lbi,fbi,gbi,iret)
985 IF(lxx.GT.0) THEN
986 IF(iret.EQ.0) THEN
987 print *,' interpolated to grid template 3.',igdtn
988 ELSEIF(iret.GT.0) THEN
989 print *,' interpolation error code ',iret
990 ENDIF
991 ENDIF
992 IF(iret.EQ.-1) iret=0
993 ENDIF
994 ENDIF
995
996 ! INVOKE MAP MASK
997 IF(iret.EQ.0) THEN
998 IF(jbk.EQ.0.AND.jb.EQ.1) THEN
999 DO i=1,mi
1000 IF(l1i(i)) THEN
1001 IF((lab.EQ.1.AND.f1i(i).LE.ab).OR. &
1002 (lab.EQ.-1.AND.f1i(i).GE.ab)) THEN
1003 ib1i=1
1004 l1i(i)=.false.
1005 ENDIF
1006 ENDIF
1007 ENDDO
1008 IF(lxx.GT.0) THEN
1009 print *,' applied post-interpolation map mask'
1010 ENDIF
1011 ELSEIF(jb.EQ.4) THEN
1012 DO i=1,mi
1013 IF(lbi(i)) THEN
1014 IF((lab.EQ.1.AND.fbi(i).LE.ab).OR. &
1015 (lab.EQ.-1.AND.fbi(i).GE.ab)) THEN
1016 ib1i=1
1017 l1i(i)=.false.
1018 ENDIF
1019 ELSE
1020 ib1i=1
1021 l1i(i)=.false.
1022 ENDIF
1023 ENDDO
1024 IF(lxx.GT.0) THEN
1025 print *,' applied fixed map mask'
1026 ENDIF
1027 ENDIF
1028
1029 ! MASK VALUES
1030 IF(lam.EQ.1.AND.ib1i.EQ.1) THEN
1031 ib1i=0
1032 DO i=1,mi
1033 IF(.NOT.l1i(i)) THEN
1034 l1i(i)=.true.
1035 f1i(i)=am
1036 IF(krv.GT.0) g1i(i)=am
1037 ENDIF
1038 ENDDO
1039 IF(lxx.GT.0) THEN
1040 print *,' substituted mask fill value ',am
1041 ENDIF
1042 ENDIF
1043 IF(lxx.GT.0) CALL instrument(3,kall3,ttot3,tmin3,tmax3)
1044
1045 ! MERGE FIELD
1046 IF(lam.EQ.5.AND.ib1i.EQ.1) THEN
1047 jdisc=gfld1%DISCIPLINE
1048 jids=-9999
1049 jjpdtn=gfld1%IPDTNUM
1050 jjpdt=-9999
1051 jjpdt(1:2)=gfld1%IPDTMPL(1:2)
1052 jjpdt(10:15)=gfld1%IPDTMPL(10:15)
1053 jgdtn=-1
1054 jgdt=-9999
1055 unpack=.true.
1056 krm=0
1057 CALL getgb2(lgm,lxm,krm,jdisc,jids,jjpdtn,jjpdt, &
1058 jgdtn,jgdt,unpack,krmx,gfldm,iret)
1059 ! CALL GETGBEM(LGM,LXM,MM,KRM,JPDS,JGDS,JENS, &
1060 ! MBUF,CBUFM,NLENM,NNUMM,MNUMM, &
1061 ! KM,KRMX,KPDSM,KGDSM,KENSM,LR,FR,IRET)
1062 km=gfldm%NGRDPTS
1063 IF(iret.EQ.0.AND.krv.GT.0) THEN
1064 gfldm%IPDTMPL(2)=gfldm%IPDTMPL(2)+1
1065 CALL getgb2(lgm,lxm,krm,gfldm%DISCIPLINE,gfldm%IDSECT, &
1066 gfldm%IPDTNUM,gfldm%IPDTMPL,gfldm%IGDTNUM, &
1067 gfldm%IGDTMPL,unpack,krmx,gfldmv,iret)
1068 ! CALL GETGBEM(LGM,LXM,MM,KRM,JPDS,JGDS,JENS, &
1069 ! MBUF,CBUFM,NLENM,NNUMM,MNUMM, &
1070 ! KM,KRMX,KPDSM,KGDSM,KENSM,LR,GR,IRET)
1071 gfldm%IPDTMPL(2)=gfldm%IPDTMPL(2)-1
1072 ENDIF
1073 IF(lxx.GT.0) THEN
1074 IF(iret.EQ.0) THEN
1075 print *,' merge field retrieved'
1076 print *,' ...PDT 4.',gfldm%IPDTNUM,'=', &
1077 (gfldm%IPDTMPL(i),i=1,gfldm%IPDTLEN)
1078 IF(krv.GT.0) &
1079 print *,' ...PDT 4.',gfldmv%IPDTNUM,'=', &
1080 (gfldmv%IPDTMPL(i),i=1,gfldmv%IPDTLEN)
1081 ELSEIF(iret.EQ.99) THEN
1082 print *,' merge field not found'
1083 ELSE
1084 print *,' merge field retrieval error code ',iret
1085 ENDIF
1086 ENDIF
1087 IF(iret.EQ.0) THEN
1088 ALLOCATE(lbi(mi),fbi(mi),gbi(mi))
1089 IF ( gfldm%IBMAP.EQ.0 .OR. gfldm%IBMAP.EQ.254 ) THEN
1090 ibm=1
1091 ELSE
1092 ibm=0
1093 ALLOCATE( gfld1%BMAP(km) ) ! dummy array
1094 ENDIF
1095 IF ( .NOT. ASSOCIATED(gfldm%LIST_OPT)) &
1096 ALLOCATE(gfldm%LIST_OPT(1))
1097 IF ( .NOT. ASSOCIATED(gfldmv%FLD)) ALLOCATE(gfldmv%FLD(km))
1098 CALL intgrib2(iv,ip,ipopt,gfldm%IGDTNUM,gfldm%IGDTMPL, &
1099 gfldm%NUM_OPT,gfldm%LIST_OPT, &
1100 km,ibm,gfldm%BMAP,gfldm%FLD,gfldmv%FLD, &
1101 igdtn,kgdti,mi, &
1102 ibbi,lbi,fbi,gbi,iret)
1103 IF(lxx.GT.0) THEN
1104 IF(iret.EQ.0) THEN
1105 print *,' interpolated to grid template 3.',igdtn
1106 ELSEIF(iret.GT.0) THEN
1107 print *,' interpolation error code ',iret
1108 ENDIF
1109 ENDIF
1110 IF(iret.EQ.-1) iret=0
1111 ENDIF
1112 IF(iret.EQ.0) THEN
1113 DO i=1,mi
1114 IF(.NOT.l1i(i).AND.lbi(i)) THEN
1115 l1i(i)=.true.
1116 f1i(i)=fbi(i)
1117 IF(krv.GT.0) g1i(i)=gbi(i)
1118 ENDIF
1119 ENDDO
1120 IF(lxx.GT.0) THEN
1121 print *,' merged output field with merge field'
1122 ENDIF
1123 ENDIF
1124 iret=0
1125 IF (ASSOCIATED(lbi)) DEALLOCATE(lbi)
1126 IF (ASSOCIATED(fbi)) DEALLOCATE(fbi)
1127 IF (ASSOCIATED(gbi)) DEALLOCATE(gbi)
1128 ENDIF
1129 IF(lxx.GT.0) CALL instrument(5,kall5,ttot5,tmin5,tmax5)
1130 ENDIF
1131
1132 ! WRITE OUTPUT FIELD
1133 IF(iret.EQ.0) THEN
1134 gfld1%IBMAP=255
1135 IF ( ib1i .EQ. 1 ) gfld1%IBMAP=0
1136 gfld1%IDRTMPL(4)=0
1137 !K5=KPDS1(5)
1138 !IDS1=KPDS1(22)
1139 !IBS1=0
1140 !NBS1=0
1141 !IF(K5.GT.0.AND.K5.LT.256) THEN
1142 ! IF(IDS(K5).GE.-128.AND.IDS(K5).LT.128) IDS1=IDS(K5)
1143 ! IF(IBS(K5).GE.-128.AND.IBS(K5).LT.128) IBS1=IBS(K5)
1144 ! IF(NBS(K5).GE.0.AND.NBS(K5).LT.256) NBS1=NBS(K5)
1145 !ENDIF
1146 !KPDS1(22)=IDS1
1147 ! Assign new GDS/GDT info to GFLD1
1148 gfld1%IGDTNUM=igdtn
1149 DEALLOCATE(gfld1%IGDTMPL)
1150 gfld1%IGDTMPL => kgdti
1151 gfld1%IGDTLEN=getgdtlen(igdtn)
1152 ! Assign new Bitmap and Data field to GFLD1
1153 gfld1%NGRDPTS=mi
1154 IF ( ASSOCIATED(gfld1%BMAP) ) DEALLOCATE(gfld1%BMAP)
1155 gfld1%BMAP => l1i
1156 IF ( ASSOCIATED(gfld1%FLD) ) DEALLOCATE(gfld1%FLD)
1157 gfld1%FLD => f1i
1158 CALL putgb2(lg2,gfld1,iret)
1159 ! CALL PUTGBEN(LG2,MI,KPDS1,KGDSI,KENS1,IBS1,NBS1,L1I,F1I,IRET)
1160 IF(iret.EQ.0) no=no+1
1161 IF(iret.EQ.0.AND.krv.GT.0) THEN
1162 IF ( ASSOCIATED(gfld1%FLD) ) DEALLOCATE(gfld1%FLD)
1163 ! Assign 2nd vector field to GFLD1
1164 gfld1%FLD => g1i
1165 ! Assign new PDT for vector field to GFLD1
1166 tmpptr => gfld1%IPDTMPL
1167 gfld1%IPDTMPL => gfldv%IPDTMPL
1168 gfld1%IDRTMPL(4)=0
1169 CALL putgb2(lg2,gfld1,iret)
1170 ! CALL PUTGBEN(LG2,MI,KPDS1,KGDSI,KENS1,IBS1,NBS1,L1I,G1I,IRET)
1171 IF(iret.EQ.0) no=no+1
1172 gfld1%IPDTMPL => tmpptr
1173 ENDIF
1174 IF(lxx.GT.0) THEN
1175 IF(iret.NE.0) THEN
1176 print *,' packing error code ',iret
1177 ELSEIF(krv.EQ.0) THEN
1178 print *,' wrote scalar field to record ',no
1179 print *,' ...PDT 4.',gfld1%IPDTNUM,'=', &
1180 (gfld1%IPDTMPL(i),i=1,gfld1%IPDTLEN)
1181 ELSE
1182 print *,' wrote vector field to records ',no-1,no
1183 print *,' ...PDT 4.',gfld1%IPDTNUM,'=', &
1184 (gfld1%IPDTMPL(i),i=1,gfld1%IPDTLEN)
1185 print *,' ...PDT 4.',gfldv%IPDTNUM,'=', &
1186 (gfldv%IPDTMPL(i),i=1,gfldv%IPDTLEN)
1187 ENDIF
1188 CALL instrument(6,kall6,ttot6,tmin6,tmax6)
1189 ENDIF
1190 ENDIF
1191
1192 ! CLEAN UP - FREE MEMORY
1193 IF (iret.EQ.0) NULLIFY(gfld1%IGDTMPL)
1194 IF (ASSOCIATED(gfld1%FLD,g1i)) then
1195 DEALLOCATE(g1i)
1196 nullify(gfld1%FLD)
1197 elseif (ALLOCATED(g1i)) then
1198 deallocate(g1i)
1199 endif
1200 CALL gf_free(gfldmv)
1201 CALL gf_free(gfldv)
1202 CALL gf_free(gfldm)
1203 CALL gf_free(gfld1)
1204 ! IF (ASSOCIATED(L1I)) DEALLOCATE(L1I)
1205 ! IF (ASSOCIATED(F1I)) DEALLOCATE(F1I)
1206 ! IF (ASSOCIATED(G1I)) DEALLOCATE(G1I)
1207
1208 END SUBROUTINE cpgb1
1209
1234 SUBROUTINE intgrib2(IV,IP,IPOPT,NGDT1,KGDT1,IDEFN1,IDEF1, &
1235 K1,IB1,L1,F1,G1, &
1236 NGDT2,KGDT2,K2,IB2,L2,F2,G2,IRET)
1237 USE gridtemplates
1238
1239 INTEGER ipopt(20)
1240 INTEGER kgdt1(*),kgdt2(*),igds(5)
1241 INTEGER idef1(idefn1)
1242 ! make idef an array to be compliant with gdt2gds. IDEFN = 0, so idef will not be used.
1243 integer idef(1)
1244 INTEGER kgds1(200),kgds2(200)
1245 LOGICAL*1 l1(k1),l2(k2)
1246 REAL f1(k1),f2(k2)
1247 REAL g1(k1),g2(k2)
1248 INTEGER kgds1f(200),kgds2f(200)
1249
1250 ! DETERMINE WHETHER INTERPOLATION IS NECESSARY
1251 IF(ip.EQ.4) THEN
1252 int=1
1253 ELSE
1254 int=0
1255 IF ( ngdt1 .EQ. ngdt2 ) THEN
1256 DO i=1,getgdtlen(ngdt1)
1257 int=max(int,abs(kgdt1(i)-kgdt2(i)))
1258 ENDDO
1259 ELSE
1260 int=1
1261 ENDIF
1262 int=min(int,1)
1263 ENDIF
1264
1265 ! COPY FIELD
1266 IF(int.EQ.0) THEN
1267 ib2=ib1
1268 DO i=1,k1
1269 l2(i)=l1(i)
1270 f2(i)=f1(i)
1271 IF(iv.NE.0) g2(i)=g1(i)
1272 ENDDO
1273 iret=-1
1274
1275 ! CONVERT GRIB2 GRID DEFINITION TEMPLATE TO GRIB1 GDS
1276 ! COMPUTE REGULARIZED GRIDS AND INTERPOLATE FIELD
1277 ELSE
1278 idefn=0
1279 igds=0
1280 igds(2)=k1
1281 igds(5)=ngdt1
1282 CALL gdt2gds(igds,kgdt1,idefn1,idef1,kgds1,igi,iret)
1283 igds(2)=k2
1284 igds(5)=ngdt2
1285 CALL gdt2gds(igds,kgdt2,idefn,idef,kgds2,igi,iret)
1286 k1f=lengdsf(kgds1,kgds1f)
1287 IF(k1f.EQ.k1) k1f=1
1288 k2f=lengdsf(kgds2,kgds2f)
1289 IF(k2f.EQ.k2) k2f=1
1290 mrl=max(k2,k2f)
1291 IF(iv.EQ.0) THEN
1292 mro=1
1293 ELSE
1294 mro=mrl
1295 ENDIF
1296 IF(k1f.GT.0.AND.k2f.GT.0) THEN
1297 CALL intgrib1(k1f,kgds1f,k2f,kgds2f,mrl,mro, &
1298 iv,ip,ipopt,kgds1,k1,ib1,l1,f1,g1,kgds2,k2, &
1299 ib2,l2,f2,g2,iret)
1300 ELSE
1301 iret=101
1302 ENDIF
1303 ENDIF
1304
1305 END SUBROUTINE intgrib2
1306
1333 SUBROUTINE intgrib1(K1F,KGDS1F,K2F,KGDS2F,MRL,MRO, &
1334 IV,IP,IPOPT,KGDS1,K1,IB1,L1,F1,G1,KGDS2,K2, &
1335 IB2,L2,F2,G2,IRET)
1336#ifdef USEIPMOD
1337 USE ip_mod
1338#endif
1339 INTEGER ipopt(20)
1340 INTEGER kgds1(200),kgds2(200)
1341 LOGICAL*1 l1(k1),l2(k2)
1342 REAL f1(k1),f2(k2),g1(k1),g2(k2)
1343 INTEGER kgds1f(200),kgds2f(200)
1344 LOGICAL*1 l1f(k1f),l2f(k2f)
1345 REAL f1f(k1f),f2f(k2f),g1f(k1f),g2f(k2f)
1346 REAL rlat(mrl),rlon(mrl),crot(mro),srot(mro)
1347
1348 ! REGLR TO REGLR SCALAR
1349 IF(k1f.EQ.1.AND.k2f.EQ.1.AND.iv.EQ.0) THEN
1350 CALL ipolates(ip,ipopt,kgds1,kgds2,k1,k2,1,ib1,l1,f1, &
1351 ki,rlat,rlon,ib2,l2,f2,iret)
1352
1353 ! IRREG TO REGLR SCALAR
1354 ELSEIF(k1f.NE.1.AND.k2f.EQ.1.AND.iv.EQ.0) THEN
1355 CALL ipxwafs2(1,k1,k1f,1, &
1356 kgds1,ib1,l1,f1,kgds1f,ib1f,l1f,f1f,iret)
1357 IF(iret.EQ.0) THEN
1358 CALL ipolates(ip,ipopt,kgds1f,kgds2,k1f,k2,1,ib1f,l1f,f1f, &
1359 ki,rlat,rlon,ib2,l2,f2,iret)
1360 ENDIF
1361
1362 ! REGLR TO IRREG SCALAR
1363 ELSEIF(k1f.EQ.1.AND.k2f.NE.1.AND.iv.EQ.0) THEN
1364 CALL ipolates(ip,ipopt,kgds1,kgds2f,k1,k2f,1,ib1,l1,f1, &
1365 ki,rlat,rlon,ib2f,l2f,f2f,iret)
1366 IF(iret.EQ.0) THEN
1367 CALL ipxwafs2(-1,k2,k2f,1, &
1368 kgds2,ib2,l2,f2,kgds2f,ib2f,l2f,f2f,iret)
1369 ENDIF
1370
1371 ! IRREG TO IRREG SCALAR
1372 ELSEIF(k1f.NE.1.AND.k2f.NE.1.AND.iv.EQ.0) THEN
1373 CALL ipxwafs2(1,k1,k1f,1, &
1374 kgds1,ib1,l1,f1,kgds1f,ib1f,l1f,f1f,iret)
1375 IF(iret.EQ.0) THEN
1376 CALL ipolates(ip,ipopt,kgds1f,kgds2f,k1f,k2f,1,ib1f,l1f,f1f, &
1377 ki,rlat,rlon,ib2f,l2f,f2f,iret)
1378 IF(iret.EQ.0) THEN
1379 CALL ipxwafs2(-1,k2,k2f,1, &
1380 kgds2,ib2,l2,f2,kgds2f,ib2f,l2f,f2f,iret)
1381 ENDIF
1382 ENDIF
1383
1384 ! REGLR TO REGLR VECTOR
1385 ELSEIF(k1f.EQ.1.AND.k2f.EQ.1.AND.iv.NE.0) THEN
1386 CALL ipolatev(ip,ipopt,kgds1,kgds2,k1,k2,1,ib1,l1,f1,g1, &
1387 ki,rlat,rlon,crot,srot,ib2,l2,f2,g2,iret)
1388 IF(iret.EQ.0.AND.ki.EQ.k2-1) THEN
1389 f2(k2)=0
1390 g2(k2)=0
1391 ENDIF
1392
1393 ! IRREG TO REGLR VECTOR
1394 ELSEIF(k1f.NE.1.AND.k2f.EQ.1.AND.iv.NE.0) THEN
1395 CALL ipxwafs2(1,k1,k1f,1, &
1396 kgds1,ib1,l1,f1,kgds1f,ib1f,l1f,f1f,iret)
1397 CALL ipxwafs2(1,k1,k1f,1, &
1398 kgds1,ib1,l1,g1,kgds1f,ib1f,l1f,g1f,iret)
1399 IF(iret.EQ.0) THEN
1400 CALL ipolatev(ip,ipopt,kgds1f,kgds2,k1f,k2,1, &
1401 ib1f,l1f,f1f,g1f, &
1402 ki,rlat,rlon,crot,srot,ib2,l2,f2,g2,iret)
1403 IF(iret.EQ.0.AND.ki.EQ.k2-1) THEN
1404 f2(k2)=0
1405 g2(k2)=0
1406 ENDIF
1407 ENDIF
1408
1409 ! REGLR TO IRREG VECTOR
1410 ELSEIF(k1f.EQ.1.AND.k2f.NE.1.AND.iv.NE.0) THEN
1411 CALL ipolatev(ip,ipopt,kgds1,kgds2f,k1,k2f,1,ib1,l1,f1,g1, &
1412 ki,rlat,rlon,crot,srot,ib2f,l2f,f2f,g2f,iret)
1413 IF(iret.EQ.0) THEN
1414 CALL ipxwafs2(-1,k2,k2f,1, &
1415 kgds2,ib2,l2,f2,kgds2f,ib2f,l2f,f2f,iret)
1416 CALL ipxwafs2(-1,k2,k2f,1, &
1417 kgds2,ib2,l2,g2,kgds2f,ib2f,l2f,g2f,iret)
1418 ENDIF
1419
1420 ! IRREG TO IRREG VECTOR
1421 ELSEIF(k1f.NE.1.AND.k2f.NE.1.AND.iv.NE.0) THEN
1422 CALL ipxwafs2(1,k1,k1f,1, &
1423 kgds1,ib1,l1,f1,kgds1f,ib1f,l1f,f1f,iret)
1424 CALL ipxwafs2(1,k1,k1f,1, &
1425 kgds1,ib1,l1,g1,kgds1f,ib1f,l1f,g1f,iret)
1426 IF(iret.EQ.0) THEN
1427 CALL ipolatev(ip,ipopt,kgds1f,kgds2f,k1f,k2f,1, &
1428 ib1f,l1f,f1f,g1f, &
1429 ki,rlat,rlon,crot,srot,ib2f,l2f,f2f,g2f,iret)
1430 IF(iret.EQ.0) THEN
1431 CALL ipxwafs2(-1,k2,k2f,1, &
1432 kgds2,ib2,l2,f2,kgds2f,ib2f,l2f,f2f,iret)
1433 CALL ipxwafs2(-1,k2,k2f,1, &
1434 kgds2,ib2,l2,g2,kgds2f,ib2f,l2f,g2f,iret)
1435 ENDIF
1436 ENDIF
1437 ENDIF
1438
1439 END SUBROUTINE intgrib1
1440
1456 FUNCTION lengdsf(KGDS,KGDSF)
1457 INTEGER kgds(200),kgdsf(200)
1458
1459 IF(kgds(1).EQ.201) THEN
1460 kgdsf=kgds
1461 lengdsf=kgds(7)*kgds(8)-kgds(8)/2
1462 ELSEIF(kgds(1).EQ.202) THEN
1463 kgdsf=kgds
1464 lengdsf=kgds(7)*kgds(8)
1465 ELSEIF(kgds(19).EQ.0.AND.kgds(20).NE.255) THEN
1466 CALL ipxwafs(1,1,1,0,kgds,dum,kgdsf,dumf,iret)
1467 IF(iret.EQ.0) THEN
1468 lengdsf=kgdsf(2)*kgdsf(3)
1469 ELSE
1470 lengdsf=0
1471 ENDIF
1472 ELSE
1473 kgdsf=kgds
1474 lengdsf=kgds(2)*kgds(3)
1475 ENDIF
1476
1477 END FUNCTION lengdsf
1478
1492 FUNCTION numpts(IGDTN,KGDT)
1493 INTEGER,INTENT(IN) :: igdtn,kgdt(*)
1494 INTEGER :: allones
1495
1496 allones = 0
1497 DO j=1,31
1498 allones=ibset(allones,j)
1499 ENDDO
1500
1501 SELECT CASE( igdtn )
1502 CASE (0:3) ! Lat/Lon
1503 IF ( kgdt(8).NE.allones .AND. kgdt(9).NE.allones ) THEN
1504 numpts = kgdt(8) * kgdt(9)
1505 ELSE
1506 numpts = -1
1507 ENDIF
1508 CASE (10) ! Mercator
1509 IF ( kgdt(8).NE.allones .AND. kgdt(9).NE.allones ) THEN
1510 numpts = kgdt(8) * kgdt(9)
1511 ELSE
1512 numpts = -1
1513 ENDIF
1514 CASE (20) ! Polar Stereographic
1515 IF ( kgdt(8).NE.allones .AND. kgdt(9).NE.allones ) THEN
1516 numpts = kgdt(8) * kgdt(9)
1517 ELSE
1518 numpts = -1
1519 ENDIF
1520 CASE (30) ! Lambert Conformal
1521 IF ( kgdt(8).NE.allones .AND. kgdt(9).NE.allones ) THEN
1522 numpts = kgdt(8) * kgdt(9)
1523 ELSE
1524 numpts = -1
1525 ENDIF
1526 CASE (40:43) ! Gaussian
1527 IF ( kgdt(8).NE.allones .AND. kgdt(9).NE.allones ) THEN
1528 numpts = kgdt(8) * kgdt(9)
1529 ELSE
1530 numpts = -1
1531 ENDIF
1532 CASE DEFAULT
1533 numpts = -1
1534 END SELECT
1535
1536 END FUNCTION numpts
1537END PROGRAM copygb2
subroutine gdt2gds(igds, igdstmpl, idefnum, ideflist, kgds, igrid, iret)
This routine converts grid information from a GRIB2 Grid Description Section as well as its Grid Defi...
Definition cnv21.F90:908
program copygb2
The command copygb2 copies all or part of one GRIB2 file to another GRIB2 file, interpolating if nece...
Definition copygb2.F90:161
subroutine cpgb(lg1, lx1, lgb, lxb, lgm, lxm, lg2, igi, kgdsi, ip, ipopt, jpds1, nuv, iuv, jpdsb, jb, jbk, lab, ab, lam, am, lxx, lwg, ids, ibs, nbs)
Copy grib files.
Definition copygb.F90:566
subroutine eusage
Print proper usage to stderr.
Definition copygb.F90:514
subroutine cpgb1(lg1, lx1, m1, cbuf1, nlen1, nnum1, mnum1, mbuf, mf, mi, igi, kgdsi, ip, ipopt, jpds1, nuv, iuv, jpdsb, jb, jbk, lab, ab, lam, am, ids, ibs, nbs, lgb, lxb, mb, cbufb, nlenb, nnumb, mnumb, lgm, lxm, mm, cbufm, nlenm, nnumm, mnumm, lg2, lxx, ks1, no, iret)
Copy one grib field.
Definition copygb.F90:787
function lengdsf(kgds, kgdsf)
Return the length of a filled grid.
Definition copygb.F90:1344
subroutine intgrib1(k1f, kgds1f, k2f, kgds2f, mrl, mro, iv, ip, ipopt, kgds1, k1, ib1, l1, f1, g1, kgds2, k2, ib2, l2, f2, g2, iret)
Interpolate field.
Definition copygb.F90:1175
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
subroutine putgb2(lugb, gfld, iret)
Pack a field into a grib2 message and write that message to a file.
Definition g2gf.F90:35
This Fortran module contains the declaration of derived type gribfield.
Definition gribmod.F90:10
This Fortran module contains info on all the available GRIB2 Grid Definition Templates used in [Secti...
integer function getgdtlen(number)
This function returns the initial length (number of entries) in the static part of specified Grid Def...