emcsfc_snow2mdl  1.5.0
 All Data Structures Files Functions Variables
snowdat.F90
Go to the documentation of this file.
1 
4 
26  module snowdat
27 
28  use program_setup, only : autosnow_file, &
29  nesdis_snow_file, &
30  nesdis_lsmask_file, &
31  afwa_snow_global_file, &
32  afwa_snow_nh_file, &
33  afwa_snow_sh_file, &
34  afwa_lsmask_nh_file, &
35  afwa_lsmask_sh_file
36 
37  use model_grid, only : imdl, &
38  jmdl
39 
40  integer :: iafwa
41  integer :: iautosnow
42  integer :: inesdis
43  integer :: jafwa
44  integer :: jautosnow
45  integer :: jnesdis
46  integer :: kgds_afwa_global(200)
48  integer :: kgds_afwa_nh(200)
50  integer :: kgds_afwa_nh_8th(200)
52  integer :: kgds_afwa_sh(200)
54  integer :: kgds_afwa_sh_8th(200)
56  integer :: kgds_autosnow(200)
57  integer :: kgds_nesdis(200)
58  integer :: mesh_nesdis
59  integer*1, allocatable :: sea_ice_nesdis(:,:)
60  logical :: bad_afwa_nh
62  logical :: bad_afwa_sh
64  logical :: bad_nesdis
66  logical :: bad_afwa_global
68  logical*1, allocatable :: bitmap_afwa_global(:,:)
70  logical*1, allocatable :: bitmap_afwa_nh(:,:)
72  logical*1, allocatable :: bitmap_afwa_sh(:,:)
74  logical*1, allocatable :: bitmap_nesdis(:,:)
75  logical*1, allocatable :: bitmap_autosnow(:,:)
77  logical :: use_nh_afwa
78  logical :: use_sh_afwa
79  logical :: use_global_afwa
80  logical :: use_autosnow
81  logical :: use_nesdis
82 
83  real :: autosnow_res
84  real :: afwa_res
85  real :: nesdis_res
86  real, allocatable :: snow_cvr_nesdis(:,:)
87  real, allocatable :: snow_cvr_autosnow(:,:)
88  real, allocatable :: snow_dep_afwa_global(:,:)
89  real, allocatable :: snow_dep_afwa_nh(:,:)
90  real, allocatable :: snow_dep_afwa_sh(:,:)
91 
92 ! the afwa 8th mesh binary data has no grib header, so set it from these
93 ! data statements. needed for ipolates routines.
94 
95  data kgds_afwa_nh_8th/5,2*512,-20826,145000,8,-80000,2*47625,0, &
96  9*0,255,180*0/
97  data kgds_afwa_sh_8th/5,2*512,20826,-125000,8,-80000,2*47625,128, &
98  9*0,255,180*0/
99  contains
117  subroutine readautosnow
118  use grib_mod ! grib 2 libraries
119 
120  implicit none
121 
122  type(gribfield) :: gfld
123 
124  integer :: iret, j, k, lugb, lugi
125  integer :: jdisc, jgdtn, jpdtn
126  integer :: jids(200), jgdt(200), jpdt(200)
127 
128  logical :: unpack
129 
130  use_autosnow = .true.
131 
132  if ( len_trim(autosnow_file) == 0 ) then
133  print*,"- WILL NOT USE AUTOSNOW DATA."
134  use_autosnow = .false.
135  return
136  end if
137 
138  print*,"- OPEN AND READ AUTOSNOW FILE ", trim(autosnow_file)
139 
140  lugb=12
141  call baopenr(lugb,autosnow_file,iret)
142 
143  if (iret /= 0) then
144  print*,'- FATAL ERROR: BAD OPEN OF FILE, IRET IS ', iret
145  call w3tage('SNOW2MDL')
146  call errexit(74)
147  endif
148 
149  call grib2_null(gfld)
150 
151  j = 0 ! search at beginning of file
152  lugi = 0 ! no grib index file
153  jdisc = 0 ! search for discipline; 0 - meteorological products
154  jpdtn = 30 ! search for product definition template number; 30 - satellite product
155  jgdtn = 0 ! search for grid definition template number; 0 - lat/lon grid
156  jids = -9999 ! array of values in identification section, set to wildcard
157  jgdt = -9999 ! array of values in grid definiation template 3.m
158  jpdt = -9999 ! array of values in product definition template 4.n
159  jpdt(1) = 1 ! search for parameter category - moisture
160  jpdt(2) = 42 ! search for parameter number - snow cover in percent.
161  unpack = .true. ! unpack data
162 
163  call getgb2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
164  unpack, k, gfld, iret)
165 
166  if (iret /=0) then
167  print*,'- FATAL ERROR: BAD DEGRIB OF FILE, IRET IS ', iret
168  call w3tage('SNOW2MDL')
169  call errexit(75)
170  endif
171 
172  print*,"- DATA VALID AT (YYYYMMDDHH): ", gfld%idsect(6),gfld%idsect(7), &
173  gfld%idsect(8),gfld%idsect(9)
174 
175  call baclose(lugb, iret)
176 
177 !-----------------------------------------------------------------------
178 ! set the grib1 kgds array from the g2 grid definition template array.
179 ! the kgds array is used by ipolates.
180 !-----------------------------------------------------------------------
181 
182  call gdt_to_gds(gfld%igdtnum, gfld%igdtmpl, gfld%igdtlen, kgds_autosnow, &
183  iautosnow, jautosnow, autosnow_res)
184 
185  allocate (bitmap_autosnow(iautosnow,jautosnow))
186  bitmap_autosnow = reshape(gfld%bmap , (/iautosnow,jautosnow/) )
187 
188  allocate (snow_cvr_autosnow(iautosnow,jautosnow))
189  snow_cvr_autosnow = reshape(gfld%fld , (/iautosnow,jautosnow/) )
190 
191  call grib2_free(gfld)
192 
193  end subroutine readautosnow
194 
222  subroutine readnesdis
223  use grib_mod
224 
225  implicit none
226 
227  integer, parameter :: iunit = 13 ! input grib file unit number
228  integer, parameter :: iunit2 = 43 ! input landmask file unit number
229 
230  integer*4, allocatable :: dummy4(:,:)
231  integer :: i, j
232  integer :: iret
233  integer :: jgds(200)
234  integer :: jpds(200)
235  integer :: lskip
236  integer, parameter :: lugi = 0 ! grib index file unit number - not used
237  integer :: jdisc, jgdtn, jpdtn, k
238  integer :: jids(200), jgdt(200), jpdt(200)
239  integer :: kgds(200)
240  integer :: kpds(200)
241  integer :: message_num
242  integer :: numbytes
243  integer :: numpts
244  integer :: isgrib
245 
246  logical :: unpack
247 
248  real, allocatable :: dummy(:,:)
249  real :: dum
250 
251  type(gribfield) :: gfld
252 
253  use_nesdis = .true.
254 
255  if ( len_trim(nesdis_snow_file) == 0 ) then
256  print*,"- WILL NOT USE NESDIS/IMS DATA."
257  use_nesdis = .false.
258  return
259  end if
260 
261  print*,"- OPEN AND READ NESDIS/IMS SNOW FILE ", trim(nesdis_snow_file)
262 
263  call grib_check(nesdis_snow_file, isgrib)
264 
265  if (isgrib==0) then
266  print*,'- FATAL ERROR: IMS FILE MUST BE GRIB 1 OR GRIB2 FORMAT'
267  call w3tage('SNOW2MDL')
268  call errexit(41)
269  end if
270 
271  call baopenr(iunit, nesdis_snow_file, iret)
272 
273  if (iret /= 0) then
274  print*,'- FATAL ERROR: BAD OPEN OF FILE, IRET IS ', iret
275  call w3tage('SNOW2MDL')
276  call errexit(73)
277  end if
278 
279  if (isgrib==1) then ! grib 1 format
280 
281 !-----------------------------------------------------------------------
282 ! tell degribber to look for requested data.
283 !-----------------------------------------------------------------------
284 
285  lskip = -1
286  jpds = -1
287  jgds = -1
288  jpds(5) = 91 ! ice cover
289  kpds = jpds
290  kgds = jgds
291 
292  print*,"- GET GRIB HEADER"
293 
294  call getgbh(iunit, lugi, lskip, jpds, jgds, numbytes, &
295  numpts, message_num, kpds, kgds, iret)
296 
297  if (iret /= 0) then
298  print*,"- FATAL ERROR: BAD DEGRIB OF HEADER. IRET IS ", iret
299  call w3tage('SNOW2MDL')
300  call errexit(72)
301  end if
302 
303  kgds_nesdis = kgds
304  inesdis = kgds(2)
305  jnesdis = kgds(3)
306 
307  mesh_nesdis = inesdis / 64
308  nesdis_res = 381. / float(mesh_nesdis) ! in km
309 
310  print*,"- DATA VALID AT (YYMMDDHH): ", kpds(8:11)
311 
312  allocate (dummy(inesdis,jnesdis))
313  allocate (sea_ice_nesdis(inesdis,jnesdis))
314  allocate (bitmap_nesdis(inesdis,jnesdis))
315 
316  print*,"- DEGRIB SEA ICE."
317 
318  call getgb(iunit, lugi, (inesdis*jnesdis), lskip, jpds, jgds, &
319  numpts, lskip, kpds, kgds, bitmap_nesdis, dummy, iret)
320 
321  if (iret /= 0) then
322  print*,"- FATAL ERROR: BAD DEGRIB OF DATA. IRET IS ", iret
323  call w3tage('SNOW2MDL')
324  call errexit(71)
325  end if
326 
327  sea_ice_nesdis = nint(dummy) ! only needed as yes/no flag
328  deallocate (dummy)
329 
330  lskip = -1
331  jpds = -1
332  jgds = -1
333  jpds(5) = 238 ! snow cover
334  kpds = jpds
335  kgds = jgds
336 
337  allocate (snow_cvr_nesdis(inesdis,jnesdis))
338 
339  print*,"- DEGRIB SNOW COVER."
340 
341  call getgb(iunit, lugi, (inesdis*jnesdis), lskip, jpds, jgds, &
342  numpts, lskip, kpds, kgds, bitmap_nesdis, snow_cvr_nesdis, iret)
343 
344  if (iret /= 0) then
345  print*,"- FATAL ERROR: BAD DEGRIB OF DATA. IRET IS ", iret
346  call w3tage('SNOW2MDL')
347  call errexit(70)
348  end if
349 
350  elseif (isgrib==2) then ! grib 2 format
351 
352  print*,"- DEGRIB SNOW COVER."
353 
354  j = 0 ! search at beginning of file
355  jdisc = 0 ! search for discipline; 0 - meteorological products
356  jpdtn = 0 ! search for product definition template number; 0 - analysis at one level
357  jgdtn = 20 ! search for grid definition template number; 20 - polar stereographic grid
358  jids = -9999 ! array of values in identification section, set to wildcard
359  jgdt = -9999 ! array of values in grid definition template 3.m
360  jpdt = -9999 ! array of values in product definition template 4.n
361  jpdt(1) = 1 ! search for parameter category - moisture
362  jpdt(2) = 201 ! search for parameter number - snow cover in percent.
363  unpack = .true. ! unpack data
364 
365  call grib2_null(gfld)
366 
367  call getgb2(iunit, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
368  unpack, k, gfld, iret)
369 
370  if (iret /=0) then
371  print*,'- FATAL ERROR: BAD DEGRIB OF FILE, IRET IS ', iret
372  call w3tage('SNOW2MDL')
373  call errexit(70)
374  endif
375 
376  print*,"- DATA VALID AT (YYYYMMDDHH): ", gfld%idsect(6),gfld%idsect(7), &
377  gfld%idsect(8),gfld%idsect(9)
378 
379 !-----------------------------------------------------------------------
380 ! set the grib1 kgds array from the g2 grid definition template array.
381 ! the kgds array is used by ipolates.
382 !-----------------------------------------------------------------------
383 
384  call gdt_to_gds(gfld%igdtnum, gfld%igdtmpl, gfld%igdtlen, kgds_nesdis, &
385  inesdis, jnesdis, dum)
386 
387  mesh_nesdis = inesdis / 64
388  nesdis_res = 381. / float(mesh_nesdis) ! in km
389 
390  if (mesh_nesdis==16) kgds_nesdis(6)=136 ! the ims 16th mesh grib2 data
391  ! is gribbed with an elliptical
392  ! earth. that is wrong. hardwire
393  ! a fix here.
394 
395  allocate (snow_cvr_nesdis(inesdis,jnesdis))
396  allocate (sea_ice_nesdis(inesdis,jnesdis))
397  allocate (bitmap_nesdis(inesdis,jnesdis))
398 
399  bitmap_nesdis = reshape(gfld%bmap , (/inesdis,jnesdis/) )
400  snow_cvr_nesdis = reshape(gfld%fld , (/inesdis,jnesdis/) )
401 
402  call grib2_free(gfld)
403 
404  print*,"- DEGRIB SEA ICE."
405 
406  j = 0 ! search at beginning of file
407  jdisc = 10 ! search for discipline; 10 - ocean products
408  jpdtn = 0 ! search for product definition template number; 0 - analysis at one level
409  jgdtn = 20 ! search for grid definition template number; 20 - polar stereographic grid
410  jids = -9999 ! array of values in identification section, set to wildcard
411  jgdt = -9999 ! array of values in grid definition template 3.m
412  jpdt = -9999 ! array of values in product definition template 4.n
413  jpdt(1) = 2 ! search for parameter category - ice
414  jpdt(2) = 0 ! search for parameter number - ice cover in percent.
415  unpack = .true. ! unpack data
416 
417  call getgb2(iunit, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
418  unpack, k, gfld, iret)
419 
420  if (iret /=0) then
421  print*,'- FATAL ERROR: BAD DEGRIB OF FILE, IRET IS ', iret
422  call w3tage('SNOW2MDL')
423  call errexit(71)
424  endif
425 
426  sea_ice_nesdis = reshape(gfld%fld , (/inesdis,jnesdis/) )
427 
428  call grib2_free(gfld)
429 
430  end if
431 
432  call baclose(iunit,iret)
433 
434 !-----------------------------------------------------------------------
435 ! the 16th mesh nesdis/ims grib data does not have a proper
436 ! bitmap section. therefore, need to read in the mask
437 ! from file. but the 96th mesh data has a proper bitmap, so use it.
438 !-----------------------------------------------------------------------
439 
440  if (mesh_nesdis == 16) then
441 
442  print*,"- OPEN NESDIS/IMS 16TH MESH LAND MASK: ", trim(nesdis_lsmask_file)
443 
444  open(iunit2, file=trim(nesdis_lsmask_file), form="formatted", &
445  iostat = iret)
446 
447  if (iret /= 0) then
448  print*,"- FATAL ERROR OPENING NESDIS/IMS LAND MASK FILE. ISTAT IS: ", iret
449  call errexit(87)
450  end if
451 
452  print*,"- READ NESDIS/IMS 16TH MESH LAND MASK."
453 
454  allocate (dummy4(inesdis,jnesdis))
455 
456  do j = 1, 1024
457  read(iunit2, 123, iostat=iret) (dummy4(i,j),i=1,1024)
458  if (iret /= 0) then
459  print*,"- FATAL ERROR READING NESDIS/IMS LAND MASK FILE. ISTAT IS: ", iret
460  call errexit(88)
461  end if
462  enddo
463 
464  close (iunit2)
465 
466 !-----------------------------------------------------------------------
467 ! the file has 0-sea, 1-land, 9-off hemi. this code expects
468 ! 0-non-land (or don't use data), 1-land (use data).
469 !-----------------------------------------------------------------------
470 
471  bitmap_nesdis=.false.
472  do j = 1, 1024
473  do i = 1, 1024
474  if (dummy4(i,j) == 1) bitmap_nesdis(i,j) = .true.
475  enddo
476  enddo
477 
478  deallocate(dummy4)
479 
480 123 FORMAT(80i1)
481 
482  endif ! is nesdis/ims data 16th mesh?
483 
484  bad_nesdis=.false.
485  call nh_climo_check(kgds_nesdis,snow_cvr_nesdis,bitmap_nesdis,inesdis,jnesdis,2,bad_nesdis)
486 
487 !-----------------------------------------------------------------------
488 ! for the 2009 nmm-b implementation, it was decided to not run with
489 ! afwa only. so even if afwa data is current and not corrupt,
490 ! but the ims is bad, then abort program. exception, if ims is very old
491 ! (there is a catastropic outage) then program will run with afwa
492 ! only. this is done by setting the nesdis_snow_file variable to
493 ! a zero length string (i.e., ims data not selected). this variable
494 ! setting is accomplished in the run script.
495 !-----------------------------------------------------------------------
496 
497  if (bad_nesdis) then
498  print*,'- FATAL ERROR: NESDIS/IMS DATA BAD, DO NOT USE.'
499  print*,'- DONT RUN PROGRAM.'
500  use_nesdis=.false.
501  call w3tage('SNOW2MDL')
502  call errexit(53)
503  stop
504  endif
505 
506  return
507 
508  end subroutine readnesdis
509 
531  subroutine readafwa
532  implicit none
533 
534  integer, parameter :: iunit=11
535  integer :: jgds(200), jpds(200), kgds(200), kpds(200)
536  integer :: istat
537  integer :: lugi, lskip, numbytes, numpts, message_num
538  integer :: isgrib
539 
540  bad_afwa_nh=.false.
541  bad_afwa_sh=.false.
542  bad_afwa_global=.false.
543 
544  use_global_afwa=.true.
545  use_nh_afwa = .true.
546  use_sh_afwa = .true.
547 
548  if (len_trim(afwa_snow_nh_file) == 0 .and. &
549  len_trim(afwa_snow_sh_file) == 0 .and. &
550  len_trim(afwa_snow_global_file) == 0) then
551  print*,"- WILL NOT USE AFWA DATA."
552  use_nh_afwa = .false.
553  use_sh_afwa = .false.
554  use_global_afwa = .false.
555  return
556  end if
557 
558  if ( len_trim(afwa_snow_global_file) > 0 ) then
559 
560  print*,"- OPEN AND READ AFWA SNOW FILE ", trim(afwa_snow_global_file)
561  call baopenr(iunit, afwa_snow_global_file, istat)
562  if (istat /= 0) then
563  print*,'- FATAL ERROR: BAD OPEN OF FILE, ISTAT IS ', istat
564  call w3tage('SNOW2MDL')
565  call errexit(60)
566  end if
567 
568 !-----------------------------------------------------------------------
569 ! tell degribber to look for requested data.
570 !-----------------------------------------------------------------------
571 
572  lugi = 0
573  lskip = -1
574  jpds = -1
575  jgds = -1
576  jpds(5) = 66 ! snow depth
577  kpds = jpds
578  kgds = jgds
579 
580  print*,"- GET GRIB HEADER"
581  call getgbh(iunit, lugi, lskip, jpds, jgds, numbytes, &
582  numpts, message_num, kpds, kgds, istat)
583 
584  if (istat /= 0) then
585  print*,"- FATAL ERROR: BAD DEGRIB OF HEADER. ISTAT IS ", istat
586  call w3tage('SNOW2MDL')
587  call errexit(61)
588  end if
589 
590  iafwa = kgds(2)
591  jafwa = kgds(3)
592  afwa_res = float(kgds(10))*0.001*111.0 ! in km.
593 
594  print*,"- DATA VALID AT (YYMMDDHH): ", kpds(8:11)
595  print*,"- DEGRIB SNOW DEPTH."
596 
597  allocate(bitmap_afwa_global(iafwa,jafwa))
598  allocate(snow_dep_afwa_global(iafwa,jafwa))
599 
600  call getgb(iunit, lugi, (iafwa*jafwa), lskip, jpds, jgds, &
601  numpts, lskip, kpds, kgds, bitmap_afwa_global, snow_dep_afwa_global, istat)
602 
603  if (istat /= 0) then
604  print*,"- FATAL ERROR: BAD DEGRIB OF DATA. ISTAT IS ", istat
605  call w3tage('SNOW2MDL')
606  call errexit(61)
607  end if
608 
609  kgds_afwa_global = kgds
610 
611  call baclose(iunit, istat)
612 
613  call nh_climo_check(kgds_afwa_global,snow_dep_afwa_global,bitmap_afwa_global,iafwa,jafwa,1,bad_afwa_global)
614 
615  if (bad_afwa_global) then
616  print*,'- WARNING: AFWA DATA BAD, DO NOT USE.'
617  use_global_afwa = .false.
618  endif
619 
620  use_nh_afwa=.false. ! use global or hemispheric files. not both.
621  use_sh_afwa=.false.
622 
623  return ! use global or hemispheric files. not both.
624 
625  else
626 
627  use_global_afwa=.false.
628 
629  endif
630 
631  if ( len_trim(afwa_snow_nh_file) > 0 ) then ! afwa nh data selected
632 
633  call grib_check(afwa_snow_nh_file, isgrib)
634 
635  if (isgrib==0) then ! old ncep binary format
636 
637  iafwa = 512
638  jafwa = 512
639  afwa_res = 47.625 ! in kilometers
640  kgds_afwa_nh = kgds_afwa_nh_8th
641 
642  allocate (snow_dep_afwa_nh(iafwa,jafwa))
643  call read_afwa_binary(afwa_snow_nh_file, snow_dep_afwa_nh)
644 
645  allocate (bitmap_afwa_nh(iafwa,jafwa))
646  call read_afwa_mask(afwa_lsmask_nh_file, bitmap_afwa_nh)
647 
648  else ! afwa data is grib
649 
650  print*,"- OPEN AND READ AFWA SNOW FILE ", trim(afwa_snow_nh_file)
651 
652  call baopenr(iunit, afwa_snow_nh_file, istat)
653 
654  if (istat /= 0) then
655  print*,'- FATAL ERROR: BAD OPEN OF FILE, ISTAT IS ', istat
656  call w3tage('SNOW2MDL')
657  call errexit(60)
658  end if
659 
660 !-----------------------------------------------------------------------
661 ! tell degribber to look for requested data.
662 !-----------------------------------------------------------------------
663 
664  lugi = 0
665  lskip = -1
666  jpds = -1
667  jgds = -1
668  jpds(5) = 66 ! snow depth
669  kpds = jpds
670  kgds = jgds
671 
672  print*,"- GET GRIB HEADER"
673  call getgbh(iunit, lugi, lskip, jpds, jgds, numbytes, &
674  numpts, message_num, kpds, kgds, istat)
675 
676  if (istat /= 0) then
677  print*,"- FATAL ERROR: BAD DEGRIB OF HEADER. ISTAT IS ", istat
678  call w3tage('SNOW2MDL')
679  call errexit(61)
680  end if
681 
682  iafwa = kgds(2)
683  jafwa = kgds(3)
684  afwa_res = float(kgds(8))*0.001 ! in km.
685 
686  print*,"- DATA VALID AT (YYMMDDHH): ", kpds(8:11)
687 
688  print*,"- DEGRIB SNOW DEPTH."
689 
690  allocate(bitmap_afwa_nh(iafwa,jafwa))
691  allocate(snow_dep_afwa_nh(iafwa,jafwa))
692 
693  call getgb(iunit, lugi, (iafwa*jafwa), lskip, jpds, jgds, &
694  numpts, lskip, kpds, kgds, bitmap_afwa_nh, snow_dep_afwa_nh, istat)
695 
696  if (istat /= 0) then
697  print*,"- FATAL ERROR: BAD DEGRIB OF DATA. ISTAT IS ", istat
698  call w3tage('SNOW2MDL')
699  call errexit(61)
700  end if
701 
702  kgds_afwa_nh = kgds
703 
704  kgds_afwa_nh(7) = -80000 ! ipolates definition of orientation angle is
705  ! 180 degrees off from grib standard.
706 
707  call baclose(iunit, istat)
708 
709  endif ! is nh afwa data grib?
710 
711  call nh_climo_check(kgds_afwa_nh,snow_dep_afwa_nh,bitmap_afwa_nh,iafwa,jafwa,1,bad_afwa_nh)
712 
713  else
714 
715  use_nh_afwa=.false.
716 
717  endif
718 
719 !-----------------------------------------------------------------------
720 ! now, read southern hemisphere data.
721 !-----------------------------------------------------------------------
722 
723  if ( len_trim(afwa_snow_sh_file) > 0 ) then
724 
725  call grib_check(afwa_snow_sh_file, isgrib)
726 
727  if (isgrib==0) then ! old ncep binary format
728 
729  iafwa = 512
730  jafwa = 512
731  afwa_res = 47.625
732  kgds_afwa_sh = kgds_afwa_sh_8th
733 
734  allocate (snow_dep_afwa_sh(iafwa,jafwa))
735  call read_afwa_binary(afwa_snow_sh_file, snow_dep_afwa_sh)
736 
737  allocate (bitmap_afwa_sh(iafwa,jafwa))
738  call read_afwa_mask(afwa_lsmask_sh_file, bitmap_afwa_sh)
739 
740  else ! sh afwa data is grib
741 
742  print*,"- OPEN AND READ AFWA SNOW FILE ", trim(afwa_snow_sh_file)
743 
744  call baopenr(iunit, afwa_snow_sh_file, istat)
745 
746  if (istat /= 0) then
747  print*,'- FATAL ERROR: BAD OPEN OF FILE, ISTAT IS ', istat
748  call w3tage('SNOW2MDL')
749  call errexit(60)
750  end if
751 
752 !-----------------------------------------------------------------------
753 ! tell degribber to look for requested data.
754 !-----------------------------------------------------------------------
755 
756  lugi = 0
757  lskip = -1
758  jpds = -1
759  jgds = -1
760  jpds(5) = 66 ! snow cover
761  kpds = jpds
762  kgds = jgds
763 
764  print*,"- GET GRIB HEADER"
765  call getgbh(iunit, lugi, lskip, jpds, jgds, numbytes, &
766  numpts, message_num, kpds, kgds, istat)
767 
768  if (istat /= 0) then
769  print*,"- FATAL ERROR: BAD DEGRIB OF HEADER. ISTAT IS ", istat
770  call w3tage('SNOW2MDL')
771  call errexit(61)
772  end if
773 
774  iafwa = kgds(2)
775  jafwa = kgds(3)
776  afwa_res = float(kgds(8))*0.001 ! in km.
777 
778  print*,"- DATA VALID AT (YYMMDDHH): ", kpds(8:11)
779 
780  print*,"- DEGRIB SNOW DEPTH."
781 
782  allocate(bitmap_afwa_sh(iafwa,jafwa))
783  allocate(snow_dep_afwa_sh(iafwa,jafwa))
784 
785  call getgb(iunit, lugi, (iafwa*jafwa), lskip, jpds, jgds, &
786  numpts, lskip, kpds, kgds, bitmap_afwa_sh, snow_dep_afwa_sh, istat)
787 
788  if (istat /= 0) then
789  print*,"- FATAL ERROR: BAD DEGRIB OF DATA. ISTAT IS ", istat
790  call w3tage('SNOW2MDL')
791  call errexit(61)
792  end if
793 
794  kgds_afwa_sh = kgds
795 
796  kgds_afwa_sh(7) = -80000 ! ipolates definition of orientation angle is
797  ! 180 degrees off from grib standard.
798 
799  call baclose(iunit, istat)
800 
801  endif ! is sh afwa data grib or not?
802 
803  call afwa_check(2)
804 
805  else
806 
807  use_sh_afwa = .false.
808 
809  endif
810 
811 !-------------------------------------------------------------------
812 !if either hemisphere is bad, don't trust all hemispheres
813 !-------------------------------------------------------------------
814 
815  if (bad_afwa_nh .or. bad_afwa_sh) then
816  print*,'- WARNING: AFWA DATA BAD, DO NOT USE.'
817  use_nh_afwa = .false.
818  use_sh_afwa = .false.
819  endif
820 
821  return
822 
823  end subroutine readafwa
824 
850  subroutine nh_climo_check(kgds_data,snow_data,bitmap_data,idata,jdata,isrc,bad)
851  use gdswzd_mod
852 
853  use program_setup, only : climo_qc_file, &
854  grib_year, grib_month, grib_day, &
855  grib_century
856 
857  use grib_mod ! for grib2 library
858 
859  implicit none
860 
861 ! describes the climo data grid.
862  integer, parameter :: iclim = 1080
863  integer, parameter :: jclim = 270
864  real, parameter :: lat11_clim = 90.0
865  real, parameter :: lon11_clim = -180.0
866  real, parameter :: dx_clim = 1./3.
867  real, parameter :: dy_clim = 1./3.
868 
869  integer, intent(in) :: idata, jdata, kgds_data(200), isrc
870  logical*1, intent(in) :: bitmap_data(idata,jdata)
871  logical, intent(out) :: bad
872  real, intent(in) :: snow_data(idata,jdata)
873 
874 ! local variables
875  integer :: idat(8), jdow, jdoy, jday
876  integer :: century, year, week, iret, lugb, i, j, ii, jj
877  integer :: lugi, jdisc, jpdtn, jgdtn, k, nret
878  integer :: jids(200), jgdt(200), jpdt(200)
879  integer :: count_nosnow_climo, count_nosnow_data
880  integer :: count_snow_climo, count_snow_data, count_grosschk_data
881 
882  logical*1, allocatable :: bitmap_clim(:,:)
883  logical :: unpack
884 
885  real, allocatable :: climo(:,:)
886  real :: fill, percent, x, y
887  real, allocatable :: xpts(:,:),ypts(:,:),rlon_data(:,:),rlat_data(:,:)
888  real :: thresh_gross, thresh
889 
890  type(gribfield) :: gfld
891 
892  bad=.false.
893  if (len_trim(climo_qc_file)==0) return
894 
895  print*,"- QC SNOW DATA IN NH."
896 
897  if (isrc==1) then
898  thresh_gross=50.0 ! afwa data is depth in meters
899  elseif (isrc==2) then
900  thresh_gross=100.0 ! nesdis/ims data is coverage in percent
901  endif
902 
903  fill=999.
904  allocate(xpts(idata,jdata))
905  allocate(ypts(idata,jdata))
906  allocate(rlon_data(idata,jdata))
907  allocate(rlat_data(idata,jdata))
908  do j=1,jdata
909  do i=1,idata
910  xpts(i,j)=i
911  ypts(i,j)=j
912  enddo
913  enddo
914 
915  print*,"- CALC LAT/LONS OF SOURCE POINTS."
916  call gdswzd(kgds_data,1,(idata*jdata),fill,xpts,ypts,rlon_data,rlat_data,nret)
917 
918  deallocate(xpts,ypts)
919 
920  if (nret /= (idata*jdata)) then
921  print*,"- WARNING: CALC FAILED. WILL NOT PERFORM QC."
922  deallocate (rlon_data,rlat_data)
923  return
924  endif
925 
926  count_grosschk_data=0
927  do j=1,jdata
928  do i=1,idata
929  if (rlat_data(i,j)>0.0 .and. bitmap_data(i,j)) then
930  if (snow_data(i,j) < 0.0 .or. snow_data(i,j) > thresh_gross) then
931  count_grosschk_data=count_grosschk_data+1
932  endif
933  endif
934  enddo
935  enddo
936 
937  if (count_grosschk_data > 1) then
938  print*,'- NUMBER OF DATA POINTS THAT FAIL GROSS CHECK ',count_grosschk_data
939  deallocate (rlon_data,rlat_data)
940  bad=.true.
941  return
942  endif
943 
944  print*,"- QC DATA SOURCE AGAINST CLIMO."
945  print*,"- OPEN CLIMO SNOW COVER FILE ",trim(climo_qc_file)
946  lugb=11
947  call baopenr(lugb,climo_qc_file,iret)
948 
949  if (iret /= 0) then
950  print*,"- WARNING: BAD OPEN, WILL NOT PERFORM QC ", iret
951  deallocate (rlon_data,rlat_data)
952  return
953  endif
954 
955 !---------------------------------------------------------------
956 ! climo file is weekly. so calculate the current week
957 ! then read that record from the climo file.
958 !---------------------------------------------------------------
959 
960  if (grib_year == 100) then
961  century = grib_century
962  else
963  century = grib_century-1
964  endif
965 
966  year = century*100 + grib_year
967 
968  idat=0
969  idat(1)=year
970  idat(2)=grib_month
971  idat(3)=grib_day
972 
973  call w3doxdat(idat,jdow,jdoy,jday)
974 
975 ! the climo file date is the beginning of the 7 day period
976 
977  week = nint((jdoy+3.)/7.)
978  if (week==0) week=52
979  if (week==53) week=1
980 
981  print*,"- READ CLIMO FOR WEEK ",week
982 
983  call grib2_null(gfld)
984 
985  j = week-1 ! search for specific week (# records to skip)
986  lugi = 0 ! no grib index file
987  jdisc = 0 ! search for discipline; 0 - meteorological products
988  jpdtn = 8 ! search for product definition template number; 8 - average
989  jgdtn = 0 ! search for grid definition template number; 0 - lat/lon grid
990  jids = -9999 ! array of values in identification section, set to wildcard
991 
992  jgdt = -9999 ! array of values in grid definition template 3.m
993  jgdt(8) = iclim ! search for assumed grid specs - i/j dimensions and corner
994  ! point lat/lons must match.
995  jgdt(9) = jclim
996  jgdt(12) = nint(lat11_clim * 1e6)
997  jgdt(13) = nint(abs(lon11_clim) * 1e6)
998 
999  jpdt = -9999 ! array of values in product definition template 4.n
1000  jpdt(1) = 1 ! search for parameter category - moisture
1001  jpdt(2) = 201 ! search for parameter number - snow cover in percent.
1002  unpack = .true. ! unpack data
1003 
1004  call getgb2(lugb, lugi, j, jdisc, jids, jpdtn, jpdt, jgdtn, jgdt, &
1005  unpack, k, gfld, iret)
1006 
1007  if (iret /= 0) then
1008  print*,"- WARNING: PROBLEM READING GRIB FILE ", iret
1009  print*,"- WILL NOT PERFORM QC."
1010  deallocate(rlon_data,rlat_data)
1011  deallocate(climo, bitmap_clim)
1012  call baclose(lugb,iret)
1013  return
1014  endif
1015 
1016  call baclose(lugb,iret)
1017 
1018  allocate(climo(iclim,jclim))
1019  climo = reshape(gfld%fld , (/iclim,jclim/) )
1020  allocate(bitmap_clim(iclim,jclim))
1021  bitmap_clim = reshape(gfld%bmap , (/iclim,jclim/) )
1022 
1023  call grib2_free(gfld)
1024 
1025 !---------------------------------------------------------------
1026 ! loop over all data points in nh. gross check data.
1027 ! afwa is a depth in meters, ims is % coverage. there should be
1028 ! no neg values or very large values. if point passes gross check,
1029 ! then check against climatology. find the
1030 ! nearest point on the climo snow cover grid. if
1031 ! climo indicates snow is likely (100% coverage), then
1032 ! check if afwa/ims has snow. if climo indicates snow is
1033 ! impossible (0% coverage), then check if afwa/ims has no snow. if
1034 ! afwa/ims differs from climo too much, then afwa/ims is
1035 ! considered suspect and will not be used.
1036 !---------------------------------------------------------------
1037 
1038  count_nosnow_climo=0
1039  count_nosnow_data=0
1040  count_snow_data=0
1041  count_snow_climo=0
1042 
1043  if (isrc==1) then
1044  thresh=.005
1045  elseif (isrc==2) then
1046  thresh=50.0
1047  endif
1048 
1049  do j=1,jdata
1050  do i=1,idata
1051  if (rlat_data(i,j)>0.0 .and. bitmap_data(i,j)) then
1052  y = (lat11_clim-rlat_data(i,j))/dy_clim + 1.0
1053  if (rlon_data(i,j)>180.0) rlon_data(i,j)=rlon_data(i,j)-360.0
1054  x = (rlon_data(i,j)-lon11_clim)/dx_clim + 1.0
1055  jj=nint(y)
1056  if (jj<1) jj=1
1057  if (jj>jclim) jj=jclim
1058  ii=nint(x)
1059  if (ii<1) ii=ii+iclim
1060  if (ii>iclim) ii=ii-iclim
1061  if (bitmap_clim(ii,jj)) then ! climo point is land
1062  if (climo(ii,jj) <1.0) then ! climo point is snow impossible
1063  count_nosnow_climo=count_nosnow_climo+1
1064  if (snow_data(i,j) == 0.0) then
1065  count_nosnow_data=count_nosnow_data+1
1066  endif
1067  endif
1068  if (climo(ii,jj) > 99.) then ! climo point is snow likely
1069  count_snow_climo=count_snow_climo+1
1070  if (snow_data(i,j) >thresh) then
1071  count_snow_data=count_snow_data+1
1072  endif
1073  endif
1074  endif
1075  endif
1076  enddo
1077  enddo
1078 
1079  percent = float(count_snow_climo-count_snow_data) / float(count_snow_climo)
1080  percent = percent*100.
1081  write(6,200) '- NUMBER OF DATA POINTS THAT SHOULD HAVE SNOW',count_snow_climo
1082  write(6,201) '- NUMBER OF THESE POINTS THAT ARE BARE GROUND',(count_snow_climo-count_snow_data), &
1083  'OR', percent, '%'
1084 
1085  200 format(1x,a45,1x,i10)
1086  201 format(1x,a45,1x,i10,1x,a2,1x,f6.2,a1)
1087 
1088  if (percent>50.0) then
1089  print*,"- WARNING: PERCENTAGE OF BARE GROUND POINTS EXCEEDS ACCEPTABLE LEVEL."
1090  print*,"- WILL NOT USE SOURCE DATA."
1091  bad=.true.
1092  endif
1093 
1094  percent = float(count_nosnow_climo-count_nosnow_data) / float(count_nosnow_climo)
1095  percent = percent*100.
1096  write(6,202) '- NUMBER OF DATA POINTS THAT SHOULD *NOT* HAVE SNOW',count_nosnow_climo
1097  write(6,203) '- NUMBER OF THESE POINTS WITH SNOW',(count_nosnow_climo-count_nosnow_data), &
1098  'OR', percent, '%'
1099 
1100  202 format(1x,a51,1x,i10)
1101  203 format(1x,a34,1x,i10,1x,a2,1x,f6.2,a1)
1102 
1103  if (percent>20.0) then
1104  print*,"- WARNING: PERCENTAGE OF POINTS WITH SNOW EXCEEDS ACCEPTABLE LEVEL."
1105  print*,"- WILL NOT USE SOURCE DATA."
1106  bad=.true.
1107  endif
1108 
1109  if (allocated(rlat_data)) deallocate (rlat_data)
1110  if (allocated(rlon_data)) deallocate (rlon_data)
1111  if (allocated(climo)) deallocate (climo)
1112  if (allocated(bitmap_clim)) deallocate (bitmap_clim)
1113 
1114  return
1115 
1116  end subroutine nh_climo_check
1117 
1122  subroutine afwa_check(hemi)
1123  use gdswzd_mod
1124 
1125  implicit none
1126 
1127  integer, intent(in) :: hemi
1128  integer :: kgds(200), nret
1129  integer, parameter :: npts=1
1130 
1131  real :: fill, xpts(npts), ypts(npts)
1132  real :: rlon(npts), rlat(npts)
1133 
1134  kgds=0
1135  fill=9999.
1136 
1137  if (hemi==1) then
1138  print*,'- QC DATA IN NH.'
1139  kgds=kgds_afwa_nh
1140  rlat=75.0
1141  rlon=-40.
1142  call gdswzd(kgds,(-1),npts,fill,xpts,ypts,rlon,rlat,nret)
1143  if (snow_dep_afwa_nh(nint(xpts(1)),nint(ypts(1))) < 0.001) then
1144  print*,'- WARNING: NO SNOW IN GREENLAND: ',snow_dep_afwa_nh(nint(xpts),nint(ypts))
1145  print*,'- DONT USE AFWA DATA.'
1146  bad_afwa_nh=.true.
1147  endif
1148  rlat=3.0
1149  rlon=-60.
1150  call gdswzd(kgds,(-1),npts,fill,xpts,ypts,rlon,rlat,nret)
1151  if (snow_dep_afwa_nh(nint(xpts(1)),nint(ypts(1))) > 0.0) then
1152  print*,'- WARNING: SNOW IN S AMERICA: ',snow_dep_afwa_nh(nint(xpts),nint(ypts))
1153  print*,'- DONT USE AFWA DATA.'
1154  bad_afwa_nh=.true.
1155  endif
1156  rlat=23.0
1157  rlon=10.
1158  call gdswzd(kgds,(-1),npts,fill,xpts,ypts,rlon,rlat,nret)
1159  if (snow_dep_afwa_nh(nint(xpts(1)),nint(ypts(1))) > 0.0) then
1160  print*,'- WARNING: SNOW IN SAHARA: ',snow_dep_afwa_nh(nint(xpts),nint(ypts))
1161  print*,'- DONT USE AFWA DATA.'
1162  bad_afwa_nh=.true.
1163  endif
1164  rlat=15.0
1165  rlon=10.
1166  call gdswzd(kgds,(-1),npts,fill,xpts,ypts,rlon,rlat,nret)
1167  if (snow_dep_afwa_nh(nint(xpts(1)),nint(ypts(1))) > 0.0) then
1168  print*,'- WARNING: SNOW IN S INDIA: ',snow_dep_afwa_nh(nint(xpts),nint(ypts))
1169  print*,'- DONT USE AFWA DATA.'
1170  bad_afwa_nh=.true.
1171  endif
1172  endif
1173 
1174  if (hemi==2) then
1175  print*,'- QC DATA IN SH.'
1176  kgds=kgds_afwa_sh
1177  rlat=-88.0
1178  rlon=0.
1179  call gdswzd(kgds,(-1),npts,fill,xpts,ypts,rlon,rlat,nret)
1180  if (snow_dep_afwa_sh(nint(xpts(1)),nint(ypts(1))) < 0.001) then
1181  print*,'- WARNING: NO SNOW IN ANTARCTICA: ',snow_dep_afwa_sh(nint(xpts),nint(ypts))
1182  print*,'- DONT USE AFWA DATA.'
1183  bad_afwa_sh=.true.
1184  endif
1185  rlat=-10.
1186  rlon=-45.
1187  call gdswzd(kgds,(-1),npts,fill,xpts,ypts,rlon,rlat,nret)
1188  if (snow_dep_afwa_sh(nint(xpts(1)),nint(ypts(1))) > 0.0) then
1189  print*,'- WARNING: SNOW IN SOUTH AMERICA: ',snow_dep_afwa_sh(nint(xpts),nint(ypts))
1190  print*,'- DONT USE AFWA DATA.'
1191  bad_afwa_sh=.true.
1192  endif
1193  rlat=-20.0
1194  rlon=130.
1195  call gdswzd(kgds,(-1),npts,fill,xpts,ypts,rlon,rlat,nret)
1196  if (snow_dep_afwa_sh(nint(xpts(1)),nint(ypts(1))) > 0.0) then
1197  print*,'- WARNING: SNOW IN AUSTRALIA: ',snow_dep_afwa_sh(nint(xpts),nint(ypts))
1198  print*,'- DONT USE AFWA DATA.'
1199  bad_afwa_sh=.true.
1200  endif
1201  rlat=-9.0
1202  rlon=25.
1203  call gdswzd(kgds,(-1),npts,fill,xpts,ypts,rlon,rlat,nret)
1204  if (snow_dep_afwa_sh(nint(xpts(1)),nint(ypts(1))) > 0.0) then
1205  print*,'- WARNING: SNOW IN AFRICA: ',snow_dep_afwa_sh(nint(xpts),nint(ypts))
1206  print*,'- DONT USE AFWA DATA.'
1207  bad_afwa_sh=.true.
1208  endif
1209  endif
1210 
1211  end subroutine afwa_check
1212 
1230  subroutine read_afwa_binary(file_name, snow_dep_afwa)
1231 
1232  implicit none
1233 
1234  character*8 :: afwa_file_info(2)
1235  character*(*), intent(in) :: file_name
1236 
1237  integer*2, allocatable :: dummy(:,:)
1238  integer :: i,j, istat
1239  integer, parameter :: iafwa = 512
1240  integer, parameter :: jafwa = 512
1241  integer, parameter :: iunit=11 ! input afwa data file
1242 
1243  real, intent(out) :: snow_dep_afwa(iafwa,jafwa)
1244 
1245  print*,"- OPEN AFWA BINARY FILE ", trim(file_name)
1246  open (iunit, file=trim(file_name), access="direct", recl=iafwa*2, iostat=istat)
1247 
1248  if (istat /= 0) then
1249  print*,'- FATAL ERROR: BAD OPEN. ISTAT IS ',istat
1250  call w3tage('SNOW2MDL')
1251  call errexit(60)
1252  end if
1253 
1254  print*,"- READ AFWA BINARY FILE ", trim(file_name)
1255  read(iunit, rec=2, iostat = istat) afwa_file_info
1256 
1257  if (istat /= 0) then
1258  print*,'- FATAL ERROR: BAD READ. ISTAT IS ',istat
1259  call w3tage('SNOW2MDL')
1260  call errexit(61)
1261  end if
1262 
1263  print*,"- AFWA DATA IS ", afwa_file_info(1), " AT TIME ", afwa_file_info(2)(2:7)
1264 
1265  allocate(dummy(iafwa,jafwa))
1266 
1267  do j = 1, jafwa
1268  read(iunit, rec=j+2, iostat=istat) (dummy(i,j),i=1,iafwa)
1269  if (istat /= 0) then
1270  print*,'- FATAL ERROR: BAD READ. ISTAT IS ',istat
1271  call w3tage('SNOW2MDL')
1272  call errexit(61)
1273  end if
1274  enddo
1275 
1276  close(iunit)
1277 
1278 !-----------------------------------------------------------------------
1279 ! "4090" is the sea ice flag. we don't use the afwa sea ice.
1280 !-----------------------------------------------------------------------
1281 
1282  where (dummy == 4090) dummy = 0
1283 
1284  snow_dep_afwa = float(dummy)
1285 
1286 !---------------------------------------------------------------------
1287 ! afwa data is a snow depth in units of tenths of inches.
1288 ! convert this to meters.
1289 !---------------------------------------------------------------------
1290 
1291  snow_dep_afwa = snow_dep_afwa * 2.54 / 1000.0
1292 
1293  deallocate (dummy)
1294 
1295  return
1296 
1297  end subroutine read_afwa_binary
1298 
1313  subroutine read_afwa_mask(file_name, bitmap_afwa)
1314  implicit none
1315 
1316  character*(*), intent(in) :: file_name
1317 
1318  integer, parameter :: iunit=11 ! input mask file
1319  integer, parameter :: iafwa = 512
1320  integer, parameter :: jafwa = 512
1321  integer :: i, j, istat
1322  integer*4, allocatable :: dummy4(:,:)
1323 
1324  logical*1, intent(out) :: bitmap_afwa(iafwa,jafwa)
1325 
1326  allocate (dummy4(iafwa,jafwa))
1327 
1328  print*,'- OPEN AFWA MASK FILE ', trim(file_name)
1329  open(iunit, file=trim(file_name), access='direct', &
1330  recl=iafwa*jafwa*4, iostat=istat)
1331 
1332  if (istat /= 0) then
1333  print*,'- FATAL ERROR: BAD OPEN. ISTAT IS ', istat
1334  call w3tage('SNOW2MDL')
1335  call errexit(62)
1336  end if
1337 
1338  print*,'- READ AFWA MASK FILE ', trim(file_name)
1339  read(iunit, rec=1, iostat=istat) dummy4
1340 
1341  if (istat /= 0) then
1342  print*,'- FATAL ERROR: BAD READ. ISTAT IS ', istat
1343  call w3tage('SNOW2MDL')
1344  call errexit(63)
1345  end if
1346 
1347  close(iunit)
1348 
1349 !-----------------------------------------------------------------------
1350 ! here -1-offhemi, 1-ocean, 2-land, 4-coast.
1351 !-----------------------------------------------------------------------
1352 
1353  bitmap_afwa = .false.
1354 
1355  do j = 1, jafwa
1356  do i = 1, iafwa
1357  if (dummy4(i,j) > 1) then
1358  bitmap_afwa(i,j) = .true.
1359  endif
1360  enddo
1361  enddo
1362 
1363  deallocate (dummy4)
1364 
1365  end subroutine read_afwa_mask
1366 
1367 
1368  end module snowdat
subroutine read_afwa_binary(file_name, snow_dep_afwa)
Read afwa binary snow depth file.
Definition: snowdat.F90:1230
subroutine grib2_null(gfld)
Nullify the grib2 gribfield pointers.
Definition: grib_utils.F90:610
subroutine readafwa
Read snow depth data and masks.
Definition: snowdat.F90:531
subroutine readnesdis
Read nesdis/ims snow cover/ice data.
Definition: snowdat.F90:222
subroutine nh_climo_check(kgds_data, snow_data, bitmap_data, idata, jdata, isrc, bad)
Check for corrupt nh snow cover data.
Definition: snowdat.F90:850
Read in data defining the model grid.
Definition: model_grid.F90:19
subroutine read_afwa_mask(file_name, bitmap_afwa)
Read afwa land mask file to get a bitmap.
Definition: snowdat.F90:1313
subroutine gdt_to_gds(igdtnum, igdstmpl, igdtlen, kgds, ni, nj, res)
Convert from the grib2 grid description template array used by the ncep grib2 library, to the grib1 grid description section array used by ncep ipolates library.
Definition: grib_utils.F90:139
subroutine grib2_free(gfld)
Deallocate the grib2 gribfield pointers.
Definition: grib_utils.F90:635
This module reads in data from the program&#39;s configuration namelist.
subroutine afwa_check(hemi)
Check for corrupt afwa data.
Definition: snowdat.F90:1122
Read and qc afwa, nesdis/ims and autosnow snow data.
Definition: snowdat.F90:26
subroutine grib_check(file_name, isgrib)
Determine whether file is grib or not.
Definition: grib_utils.F90:23
subroutine readautosnow
Read autosnow snow cover.
Definition: snowdat.F90:117