emcsfc_snow2mdl  1.9.0
snow2mdl.F90
Go to the documentation of this file.
1 
4 
20  module snow2mdl
21 
22  use program_setup, only : lat_threshold, &
26  grib_day, &
27  grib_century, &
28  grib_hour, &
29  grib_month, &
30  grib_year, &
32 
33  use model_grid, only : resol_mdl, &
34  imdl, &
35  jmdl, &
36  ijmdl, &
37  ipts_mdl, &
38  jpts_mdl, &
39  lsmask_mdl, &
41  lats_mdl, &
42  lat11, latlast, &
43  lon11, lonlast, &
44  lons_mdl, &
45  kgds_mdl, &
46  grid_id_mdl, &
47  thinned, &
49 
50  use snowdat, only : nesdis_res, &
51  afwa_res, &
52  autosnow_res, &
53  inesdis, &
54  jnesdis, &
55  mesh_nesdis, &
58  bitmap_nesdis, &
59  iafwa, &
60  jafwa, &
67  iautosnow, &
68  jautosnow, &
72  use_nh_afwa, &
73  use_sh_afwa, &
74  use_nesdis, &
75  use_autosnow, &
76  kgds_nesdis, &
78  kgds_afwa_nh, &
79  kgds_afwa_sh, &
80  kgds_autosnow, &
82 
83  private
84 
85  real, allocatable :: snow_cvr_mdl(:,:)
86  real, allocatable :: snow_dep_mdl(:,:)
87 
88  public :: interp
89 
90  contains
157  subroutine interp
158  use gdswzd_mod
159 
160  implicit none
161 
162  integer :: i, j, ii, jj, ij
163  integer :: ijmdl2, istart, iend, imid, iii
164  integer, allocatable :: idum(:,:)
165  integer :: int_opt, ipopt(20)
166  integer :: kgds_mdl_tmp(200)
167  integer :: no, ibo, iret, nret
168 
169  logical*1, allocatable :: bitmap_mdl(:)
170 
171  real :: gridi(1)
172  real :: gridj(1)
173  real, allocatable :: lsmask_1d(:)
174  real, allocatable :: snow_cvr_mdl_1d(:)
175  real, allocatable :: snow_dep_mdl_tmp(:)
176  real :: sumc, sumd, x1, r, fraction, gridis, gridie
177  real, parameter :: undefined_value = -999.
178 
179 !----------------------------------------------------------------------
180 ! for model grids fully or partially located in the southern
181 ! hemisphere, the user must select sh data (either autosnow or afwa).
182 ! this restriction is relaxed for the ndas domain, which is
183 ! has a small part below the equator.
184 !----------------------------------------------------------------------
185 
186  if (minval(lats_mdl) < -10.0 .and. .not. use_global_afwa .and. .not. use_sh_afwa .and. .not. use_autosnow) then
187  print*,"- FATAL ERROR: MUST SELECT EITHER AFWA OR AUTOSNOW DATA FOR MODEL GRID WITH SH POINTS."
188  call w3tage('SNOW2MDL')
189  call errexit(54)
190  endif
191 
192 !----------------------------------------------------------------------
193 ! if model grid is totally within the southern hemisphere, set flags
194 ! so that nh afwa and nh nesdis/ims is not processed. these flags
195 ! are set to false if user inadvertantly selects these data sources.
196 !----------------------------------------------------------------------
197 
198  if (maxval(lats_mdl) < 0.0) then
199  use_nh_afwa=.false.
200  use_nesdis=.false.
201  endif
202 
203 !----------------------------------------------------------------------
204 ! if model grid is partially or fully located in the northern
205 ! hemisphere, the user must select either nesdis/ims or nh afwa data.
206 !----------------------------------------------------------------------
207 
208  if (maxval(lats_mdl) > 0.0 .and. .not. use_global_afwa .and. .not. use_nh_afwa .and. .not. use_nesdis) then
209  print*,"- FATAL ERROR: MUST SELECT EITHER NESDIS/IMS OR AFWA DATA FOR MODEL GRID WITH NH POINTS."
210  call w3tage('SNOW2MDL')
211  call errexit(54)
212  endif
213 
214 !----------------------------------------------------------------------
215 ! if model grid is totally within the northern hemisphere, set flag
216 ! so that sh data is not processed. these flags are set to false
217 ! if user inadvertantly selects this data.
218 !----------------------------------------------------------------------
219 
220  if (minval(lats_mdl) > 0.0) then ! is model grid totally within northern hemisphere?
221  use_sh_afwa=.false.
222  use_autosnow=.false.
223  endif
224 
225 !----------------------------------------------------------------------
226 ! if selected, interpolate nesdis/ims data to model grid.
227 !----------------------------------------------------------------------
228 
229  nesdis_ims : if (use_nesdis) then
230 
231  ipopt = 0
232  if (nesdis_res < (0.5*resol_mdl)) then
233  print*,"- INTERPOLATE NH NESDIS/IMS DATA TO MODEL GRID USING BUDGET METHOD."
234  ipopt(1)=2 ! break model grid cell into 25 points.
235  ipopt(2:4)=1 ! 25 points are weighted equally.
236  ipopt(5)=10 ! 10% coverage of valid data in box
237  ipopt(20) = nint(100.0 / nesdis_res) + 1 ! search box width of 100 km.
238  kgds_mdl_tmp = kgds_mdl
239  kgds_mdl_tmp(1) = kgds_mdl_tmp(1) - 255 ! subset of grid
240  int_opt = 3
241  no = ijmdl
242  else
243  print*,"- INTERPOLATE NH NESDIS/IMS DATA TO MODEL GRID USING NEIGHBOR METHOD."
244  ipopt(1) = nint(100.0 / nesdis_res) + 1 ! search box width of 100 km.
245  kgds_mdl_tmp = kgds_mdl
246  kgds_mdl_tmp(1) = -1 ! for subsection of model grid.
247  int_opt = 2
248  no = ijmdl ! an input when kgds(1) < 0 (subset of grid)
249  end if
250 
251  allocate (snow_cvr_mdl_1d(ijmdl))
252  snow_cvr_mdl_1d = 0.0
253 
254  allocate (bitmap_mdl(ijmdl))
255  bitmap_mdl=.false. ! if interpolation routine can't find data
256  ! at a point, this flag is false.
257 
258  call ipolates(int_opt, ipopt, kgds_nesdis, kgds_mdl_tmp, &
259  (inesdis*jnesdis), ijmdl, &
261  no, lats_mdl, lons_mdl, ibo, bitmap_mdl, &
262  snow_cvr_mdl_1d, iret)
263 
264  deallocate (bitmap_nesdis, snow_cvr_nesdis)
265 
266  if (iret /= 0) then
267  print*,"- FATAL ERROR: IN INTERPOLATION ROUTINE. IRET IS: ", iret
268  call w3tage('SNOW2MDL')
269  call errexit(55)
270  endif
271 
272 !-----------------------------------------------------------------------
273 ! if the interpolation routines did not find valid nesdis/ims data
274 ! in the vicinity of the model point, need to set a default value
275 ! of snow cover. south of user-defined latitude threshold, set
276 ! to zero. otherwise, see if the nearest neighbor nesdis/ims point is
277 ! sea ice. if so, assume model point is snow covered.
278 !-----------------------------------------------------------------------
279 
280  do ij = 1, ijmdl
281  if (lats_mdl(ij) < 0.0) cycle ! only consider nh model points
282  if (.not. bitmap_mdl(ij)) then
283  if (lats_mdl(ij) <= lat_threshold) then
284  snow_cvr_mdl_1d(ij) = 0.0
285  else
286  call gdswzd(kgds_nesdis,-1,1,undefined_value,gridi,gridj, &
287  lons_mdl(ij),lats_mdl(ij),nret)
288  if (nret /= 1) then
289  print*,"- WARNING: MODEL POINT OUTSIDE NESDIS/IMS GRID: ", ipts_mdl(ij), jpts_mdl(ij)
290  snow_cvr_mdl_1d(ij) = 0.0
291  else
292  ii = nint(gridi(1))
293  jj = nint(gridj(1))
294  if (sea_ice_nesdis(ii,jj) == 1) then
295  snow_cvr_mdl_1d(ij) = 100.0
296  else
297  snow_cvr_mdl_1d(ij) = 0.0
298  end if
299  end if
300  end if
301  end if
302  enddo
303 
304  deallocate (sea_ice_nesdis)
305  deallocate (bitmap_mdl)
306 
307  endif nesdis_ims
308 
309 !----------------------------------------------------------------------
310 ! now interpolate nh afwa snow depth data.
311 !----------------------------------------------------------------------
312 
313  global_afwa : if (use_global_afwa) then
314 
315 !----------------------------------------------------------------------
316 ! determine interpolation method based on the resolution of
317 ! afwa data and the model grid.
318 !----------------------------------------------------------------------
319 
320  ipopt = 0
321  if (afwa_res < (0.5*resol_mdl)) then
322  print*,"- INTERPOLATE GLOBAL AFWA DATA TO MODEL GRID USING BUDGET METHOD."
323  ipopt(1)=-1 ! break model grid cell into 25 points.
324  ipopt(2)=-1 ! 25 points are weighted equally.
325  ipopt(20) = nint(100.0 / afwa_res) + 1 ! search box width of 100 km.
326  kgds_mdl_tmp = kgds_mdl
327  kgds_mdl_tmp(1) = kgds_mdl_tmp(1) - 255 ! subset of grid
328  no = ijmdl
329  int_opt = 3
330  else
331  print*,"- INTERPOLATE GLOBAL AFWA DATA TO MODEL GRID USING NEIGHBOR METHOD."
332  ipopt(1) = nint(100.0 / afwa_res) + 1 ! search box width of 100 km.
333  kgds_mdl_tmp = kgds_mdl
334  kgds_mdl_tmp(1) = -1 ! for subsection of model grid.
335  int_opt = 2
336  no = ijmdl ! an input when kgds(1) < 0 (subset of grid)
337  end if
338 
339  allocate (snow_dep_mdl_tmp(ijmdl))
340  snow_dep_mdl_tmp = 0.0
341 
342  allocate (bitmap_mdl(ijmdl))
343  bitmap_mdl = .false.
344 
345  call ipolates(int_opt, ipopt, kgds_afwa_global, kgds_mdl_tmp, &
346  (iafwa*jafwa), ijmdl, &
348  no, lats_mdl, lons_mdl, ibo, bitmap_mdl, &
349  snow_dep_mdl_tmp, iret)
350 
352 
353  if (iret /= 0) then
354  print*,"- FATAL ERROR IN INTERPOLATION ROUTINE. IRET IS: ", iret
355  call w3tage('SNOW2MDL')
356  call errexit(55)
357  endif
358 
359 !----------------------------------------------------------------------
360 ! if interpolation did not find afwa data near the model point, then
361 ! use a nominal value based on latitude threshold. this value
362 ! may be overwritten below depending on nesdis/ims cover data (if
363 ! user selects nesdis/ims data).
364 !----------------------------------------------------------------------
365 
366  do ij = 1, ijmdl
367  if (.not. bitmap_mdl(ij)) then
368  if (abs(lats_mdl(ij)) >= lat_threshold) then
369  snow_dep_mdl_tmp(ij) = min_snow_depth
370  else
371  snow_dep_mdl_tmp(ij) = 0.0
372  endif
373  endif
374  enddo
375 
376  deallocate(bitmap_mdl)
377 
378  endif global_afwa
379 
380 !----------------------------------------------------------------------
381 ! now interpolate nh afwa snow depth data.
382 !----------------------------------------------------------------------
383 
384  nh_afwa : if (use_nh_afwa) then
385 
386 !----------------------------------------------------------------------
387 ! determine interpolation method based on the resolution of
388 ! afwa data and the model grid.
389 !----------------------------------------------------------------------
390 
391  ipopt = 0
392  if (afwa_res < (0.5*resol_mdl)) then
393  print*,"- INTERPOLATE NH AFWA DATA TO MODEL GRID USING BUDGET METHOD."
394  ipopt(1)=-1 ! break model grid cell into 25 points.
395  ipopt(2)=-1 ! 25 points are weighted equally.
396  ipopt(20) = nint(100.0 / afwa_res) + 1 ! search box width of 100 km.
397  kgds_mdl_tmp = kgds_mdl
398  kgds_mdl_tmp(1) = kgds_mdl_tmp(1) - 255 ! subset of grid
399  no = ijmdl
400  int_opt = 3
401  else
402  print*,"- INTERPOLATE NH AFWA DATA TO MODEL GRID USING NEIGHBOR METHOD."
403  ipopt(1) = nint(100.0 / afwa_res) + 1 ! search box width of 100 km.
404  kgds_mdl_tmp = kgds_mdl
405  kgds_mdl_tmp(1) = -1 ! for subsection of model grid.
406  int_opt = 2
407  no = ijmdl ! an input when kgds(1) < 0 (subset of grid)
408  end if
409 
410  allocate (snow_dep_mdl_tmp(ijmdl))
411  snow_dep_mdl_tmp = 0.0
412 
413  allocate (bitmap_mdl(ijmdl))
414  bitmap_mdl = .false.
415 
416  call ipolates(int_opt, ipopt, kgds_afwa_nh, kgds_mdl_tmp, &
417  (iafwa*jafwa), ijmdl, &
419  no, lats_mdl, lons_mdl, ibo, bitmap_mdl, &
420  snow_dep_mdl_tmp, iret)
421 
422  deallocate(bitmap_afwa_nh, snow_dep_afwa_nh)
423 
424  if (iret /= 0) then
425  print*,"- FATAL ERROR IN INTERPOLATION ROUTINE. IRET IS: ", iret
426  call w3tage('SNOW2MDL')
427  call errexit(55)
428  endif
429 
430 !----------------------------------------------------------------------
431 ! if interpolation did not find afwa data near the model point, then
432 ! use a nominal value based on latitude threshold. this value
433 ! may be overwritten below depending on nesdis/ims cover data (if
434 ! user selects nesdis/ims data).
435 !----------------------------------------------------------------------
436 
437  do ij = 1, ijmdl
438  if (lats_mdl(ij) >= 0.) then ! only consider model pts in n hemi.
439  if (.not. bitmap_mdl(ij)) then
440  if (abs(lats_mdl(ij)) >= lat_threshold) then
441  snow_dep_mdl_tmp(ij) = min_snow_depth
442  else
443  snow_dep_mdl_tmp(ij) = 0.0
444  endif
445  endif
446  endif
447  enddo
448 
449  deallocate(bitmap_mdl)
450 
451  endif nh_afwa
452 
453 !----------------------------------------------------------------------
454 ! if nh data selected, use it to determine the cover and depth
455 ! on the model grid.
456 !----------------------------------------------------------------------
457 
458  allocate (snow_dep_mdl(imdl,jmdl))
459  allocate (snow_cvr_mdl(imdl,jmdl))
460  snow_cvr_mdl = 0.0
461  snow_dep_mdl = 0.0
462 
463  if (use_global_afwa .and. use_nesdis) then ! set depth/cover on nesdis ims/afwa blend
464  print*,"- BLEND NESDIS/IMS AND AFWA DATA IN NH."
465  do ij = 1, ijmdl
466  if (lats_mdl(ij) >= 0.0) then
467  if (snow_cvr_mdl_1d(ij) >= snow_cvr_threshold) then
468  snow_dep_mdl(ipts_mdl(ij),jpts_mdl(ij)) = &
469  max(snow_dep_mdl_tmp(ij), min_snow_depth)
470  endif
471  snow_cvr_mdl(ipts_mdl(ij),jpts_mdl(ij)) = snow_cvr_mdl_1d(ij)
472  endif
473  enddo
474  deallocate (snow_cvr_mdl_1d)
475  elseif (use_nh_afwa .and. use_nesdis) then ! set depth/cover on nesdis ims/afwa blend
476  print*,"- BLEND NESDIS/IMS AND AFWA DATA IN NH."
477  do ij = 1, ijmdl
478  if (lats_mdl(ij) >= 0.0) then
479  if (snow_cvr_mdl_1d(ij) >= snow_cvr_threshold) then
480  snow_dep_mdl(ipts_mdl(ij),jpts_mdl(ij)) = &
481  max(snow_dep_mdl_tmp(ij), min_snow_depth)
482  endif
483  snow_cvr_mdl(ipts_mdl(ij),jpts_mdl(ij)) = snow_cvr_mdl_1d(ij)
484  endif
485  enddo
486  deallocate (snow_cvr_mdl_1d)
487  deallocate (snow_dep_mdl_tmp)
488  elseif (use_global_afwa) then ! set depth/cover on afwa only
489  print*,"- SET DEPTH/COVER FROM AFWA DATA IN NH."
490  do ij = 1, ijmdl
491  if (lats_mdl(ij) >= 0.0) then
492  if (snow_dep_mdl_tmp(ij) > 0.0) then
493  snow_cvr_mdl(ipts_mdl(ij),jpts_mdl(ij)) = 100.0
494  snow_dep_mdl(ipts_mdl(ij),jpts_mdl(ij)) = snow_dep_mdl_tmp(ij)
495  endif
496  endif
497  enddo
498  elseif (use_nh_afwa) then ! set depth/cover on afwa only
499  print*,"- SET DEPTH/COVER FROM AFWA DATA IN NH."
500  do ij = 1, ijmdl
501  if (lats_mdl(ij) >= 0.0) then
502  if (snow_dep_mdl_tmp(ij) > 0.0) then
503  snow_cvr_mdl(ipts_mdl(ij),jpts_mdl(ij)) = 100.0
504  snow_dep_mdl(ipts_mdl(ij),jpts_mdl(ij)) = snow_dep_mdl_tmp(ij)
505  endif
506  endif
507  enddo
508  deallocate (snow_dep_mdl_tmp)
509  elseif (use_nesdis) then ! set depth/cover on nesdis/ims only
510  print*,"- SET DEPTH/COVER FROM NESDIS/IMS DATA IN NH."
511  do ij = 1, ijmdl
512  if (lats_mdl(ij) >= 0.0) then
513  if (snow_cvr_mdl_1d(ij) >= snow_cvr_threshold) then
515  endif
516  snow_cvr_mdl(ipts_mdl(ij),jpts_mdl(ij)) = snow_cvr_mdl_1d(ij)
517  endif
518  enddo
519  deallocate (snow_cvr_mdl_1d)
520  end if
521 
522 !----------------------------------------------------------------------
523 ! if selected, interpolate autosnow data to model grid.
524 !----------------------------------------------------------------------
525 
526  autosnow : if (use_autosnow) then
527 
528  ipopt = 0
529  if (autosnow_res < (0.5*resol_mdl)) then
530  print*,"- INTERPOLATE AUTOSNOW DATA TO MODEL GRID USING BUDGET METHOD."
531  ipopt(1)=2 ! break model grid cell into 25 points.
532  ipopt(2:4)=1 ! 25 points are weighted equally.
533  ipopt(5)=10 ! 10% coverage of valid data in box
534  ipopt(20) = nint(100.0 / autosnow_res) + 1 ! search box width of 100 km.
535  kgds_mdl_tmp = kgds_mdl
536  kgds_mdl_tmp(1) = kgds_mdl_tmp(1) - 255 ! subset of grid
537  int_opt = 3
538  no = ijmdl
539  else
540  print*,"- INTERPOLATE AUTOSNOW DATA TO MODEL GRID USING NEIGHBOR METHOD."
541  ipopt(1) = nint(100.0 / autosnow_res) + 1 ! search box width of 100 km.
542  kgds_mdl_tmp = kgds_mdl
543  kgds_mdl_tmp(1) = -1 ! for subsection of model grid.
544  int_opt = 2
545  no = ijmdl ! an input when kgds(1) < 0 (subset of grid)
546  end if
547 
548  allocate (snow_cvr_mdl_1d(ijmdl))
549  snow_cvr_mdl_1d = 0.0
550 
551  allocate (bitmap_mdl(ijmdl))
552  bitmap_mdl=.false. ! if interpolation routine can't find data
553  ! at a point, this flag is false.
554 
555  call ipolates(int_opt, ipopt, kgds_autosnow, kgds_mdl_tmp, &
558  no, lats_mdl, lons_mdl, ibo, bitmap_mdl, &
559  snow_cvr_mdl_1d, iret)
560 
561  deallocate (snow_cvr_autosnow, bitmap_autosnow)
562 
563  if (iret /= 0) then
564  print*,"- FATAL ERROR IN INTERPOLATION ROUTINE. IRET IS: ", iret
565  call w3tage('SNOW2MDL')
566  call errexit(55)
567  endif
568 
569 !----------------------------------------------------------------------
570 ! if interpolation fails to find autosnow data, set the cover
571 ! at the model point to a nomimal value.
572 !----------------------------------------------------------------------
573 
574  do ij = 1, ijmdl
575  if (lats_mdl(ij) < 0.0) then
576  if (.not. bitmap_mdl(ij)) then
577  if (abs(lats_mdl(ij)) <= lat_threshold) then
578  snow_cvr_mdl_1d(ij) = 0.0
579  else
580  snow_cvr_mdl_1d(ij) = 100.0
581  end if
582  end if
583  end if
584  enddo
585 
586  deallocate (bitmap_mdl)
587 
588  endif autosnow
589 
590 !----------------------------------------------------------------------
591 ! now interpolate sh afwa snow depth data.
592 !----------------------------------------------------------------------
593 
594  sh_afwa : if (use_sh_afwa) then
595 
596 !----------------------------------------------------------------------
597 ! determine interpolation method based on the resolution of
598 ! afwa data and the model grid.
599 !----------------------------------------------------------------------
600 
601  ipopt = 0
602  if (afwa_res < (0.5*resol_mdl)) then
603  print*,"- INTERPOLATE SH AFWA DATA TO MODEL GRID USING BUDGET METHOD."
604  ipopt(1)=-1 ! break model grid cell into 25 points.
605  ipopt(2)=-1 ! 25 points are weighted equally.
606  ipopt(20) = nint(100.0 / afwa_res) + 1 ! search box width of 100 km.
607  kgds_mdl_tmp = kgds_mdl
608  kgds_mdl_tmp(1) = kgds_mdl_tmp(1) - 255 ! subset of grid
609  no = ijmdl
610  int_opt = 3
611  else
612  print*,"- INTERPOLATE SH AFWA DATA TO MODEL GRID USING NEIGHBOR METHOD."
613  ipopt(1) = nint(100.0 / afwa_res) + 1 ! search box width of 100 km.
614  kgds_mdl_tmp = kgds_mdl
615  kgds_mdl_tmp(1) = -1 ! for subsection of model grid.
616  int_opt = 2
617  no = ijmdl ! an input when kgds(1) < 0 (subset of grid)
618  end if
619 
620  allocate (snow_dep_mdl_tmp(ijmdl))
621  snow_dep_mdl_tmp = 0.0
622 
623  allocate (bitmap_mdl(ijmdl))
624  bitmap_mdl = .false.
625 
626  call ipolates(int_opt, ipopt, kgds_afwa_sh, kgds_mdl_tmp, &
627  (iafwa*jafwa), ijmdl, &
629  no, lats_mdl, lons_mdl, ibo, bitmap_mdl, &
630  snow_dep_mdl_tmp, iret)
631 
632  if (iret /= 0) then
633  print*,"- FATAL ERROR IN INTERPOLATION ROUTINE. IRET IS: ", iret
634  call w3tage('SNOW2MDL')
635  call errexit(55)
636  endif
637 
638  deallocate (bitmap_afwa_sh, snow_dep_afwa_sh)
639 
640 !----------------------------------------------------------------------
641 ! if interpolation does not find afwa data, set model point to
642 ! a nominal value.
643 !----------------------------------------------------------------------
644 
645  do ij = 1, ijmdl
646  if (lats_mdl(ij) < 0.) then
647  if (.not. bitmap_mdl(ij)) then
648  if (abs(lats_mdl(ij)) >= lat_threshold) then
649  snow_dep_mdl_tmp(ij) = min_snow_depth
650  else
651  snow_dep_mdl_tmp(ij) = 0.0
652  endif
653  endif
654  endif
655  enddo
656 
657  deallocate(bitmap_mdl)
658 
659  endif sh_afwa
660 
661 !----------------------------------------------------------------------
662 ! if sh data selected, use it to determine the cover and depth
663 ! on the model grid.
664 !----------------------------------------------------------------------
665 
666  if ((use_sh_afwa .or. use_global_afwa) .and. use_autosnow) then ! set depth/cover on autosnow/afwa blend
667  print*,"- BLEND AUTOSNOW AND AFWA DATA IN SH."
668  do ij = 1, ijmdl
669  if (lats_mdl(ij) < 0.0) then
670  if (snow_cvr_mdl_1d(ij) >= snow_cvr_threshold) then
671  snow_dep_mdl(ipts_mdl(ij),jpts_mdl(ij)) = &
672  max(snow_dep_mdl_tmp(ij), min_snow_depth)
673  endif
674  snow_cvr_mdl(ipts_mdl(ij),jpts_mdl(ij)) = snow_cvr_mdl_1d(ij)
675  endif
676  enddo
677  deallocate (snow_cvr_mdl_1d)
678  deallocate (snow_dep_mdl_tmp)
679  elseif (use_sh_afwa .or. use_global_afwa) then ! set depth/cover on afwa only
680  print*,"- SET DEPTH/COVER FROM AFWA DATA IN SH."
681  do ij = 1, ijmdl
682  if (lats_mdl(ij) < 0.0) then
683  if (snow_dep_mdl_tmp(ij) > 0.0) then
684  snow_cvr_mdl(ipts_mdl(ij),jpts_mdl(ij)) = 100.0
685  snow_dep_mdl(ipts_mdl(ij),jpts_mdl(ij)) = snow_dep_mdl_tmp(ij)
686  endif
687  endif
688  enddo
689  deallocate (snow_dep_mdl_tmp)
690  elseif (use_autosnow) then ! set depth/cover on autosnow only
691  print*,"- SET DEPTH/COVER FROM AUTOSNOW IN SH."
692  do ij = 1, ijmdl
693  if (lats_mdl(ij) < 0.0) then
694  if (snow_cvr_mdl_1d(ij) >= snow_cvr_threshold) then
696  endif
697  snow_cvr_mdl(ipts_mdl(ij),jpts_mdl(ij)) = snow_cvr_mdl_1d(ij)
698  endif
699  enddo
700  deallocate (snow_cvr_mdl_1d)
701  end if
702 
703 !----------------------------------------------------------------------
704 ! if a global model grid, and if running on thinned grid, then
705 ! take a linear weighting of full points located within the thin points.
706 ! "4" is grid indicator for a gaussian grid.
707 !----------------------------------------------------------------------
708 
709  if (kgds_mdl(1) == 4 .and. thinned) then
710 
711  ijmdl2 = sum(lonsperlat_mdl) * 2
712  allocate (snow_cvr_mdl_1d(ijmdl2))
713  allocate (lsmask_1d(ijmdl2))
714  allocate (snow_dep_mdl_tmp(ijmdl2))
715 
716  lsmask_1d = 0.0
717  snow_cvr_mdl_1d = 0.0
718  snow_dep_mdl_tmp = 0.0
719 
720  ij = 0
721  do j = 1, jmdl
722  jj = j
723  if (jj > jmdl/2) jj = jmdl - j + 1
724  r = float(imdl) / float(lonsperlat_mdl(jj))
725  do i = 1, lonsperlat_mdl(jj)
726  ij = ij + 1
727  x1 = (i-1)*r
728  imid = nint(x1+1.0)
729  lsmask_1d(ij) = lsmask_mdl_sav(imid,j)
730  if (lsmask_mdl_sav(imid,j) == 0.0) cycle
731  gridis=x1+1.0-r/2.
732  istart=nint(gridis)
733  gridie=x1+1.0+r/2.
734  iend=nint(gridie)
735  sumc = 0.0 ! cover
736  sumd = 0.0 ! depth
737  do ii = istart, iend
738  if (ii == istart) then
739  fraction = 0.5 - (gridis - float(istart))
740  elseif (ii == iend) then
741  fraction = 0.5 + (gridie - float(iend))
742  else
743  fraction = 1.0
744  endif
745  if (fraction < 0.0001) cycle
746  iii = ii
747  if (iii < 1) iii = imdl + iii
748  sumc = sumc + fraction * snow_cvr_mdl(iii,j)
749  sumd = sumd + fraction * snow_dep_mdl(iii,j)
750  enddo
751  snow_cvr_mdl_1d(ij) = sumc / r
752  snow_dep_mdl_tmp(ij) = 0.0
753  if (snow_cvr_mdl_1d(ij) > snow_cvr_threshold) then
754  snow_dep_mdl_tmp(ij) = max(sumd / r,min_snow_depth)
755  end if
756  enddo
757  enddo
758 
759  deallocate (lsmask_mdl_sav)
760 
761 !----------------------------------------------------------------------
762 ! now place thinned points into 2-d array for output.
763 !----------------------------------------------------------------------
764 
765  allocate (idum(imdl,jmdl))
766  idum = 0
767  call uninterpred(1, idum, lsmask_1d, lsmask_mdl, imdl, jmdl, ijmdl2, lonsperlat_mdl)
768  call uninterpred(1, idum, snow_cvr_mdl_1d, snow_cvr_mdl, imdl, jmdl, ijmdl2, lonsperlat_mdl)
769  deallocate(snow_cvr_mdl_1d)
770  call uninterpred(1, idum, snow_dep_mdl_tmp, snow_dep_mdl, imdl, jmdl, ijmdl2, lonsperlat_mdl)
771  deallocate(snow_dep_mdl_tmp)
772  deallocate(idum)
773 
774  end if
775 
776 !----------------------------------------------------------------------
777 ! grib the interpolated data.
778 !----------------------------------------------------------------------
779 
780  if (output_grib2) then
781  print*,"- OUTPUT SNOW ANALYSIS DATA IN GRIB2 FORMAT"
782  call write_grib2
783  else
784  print*,"- OUTPUT SNOW ANALYSIS DATA IN GRIB1 FORMAT"
785  call write_grib1
786  endif
787 
788  deallocate (snow_cvr_mdl)
789  deallocate (snow_dep_mdl)
790 
791  return
792 
793  end subroutine interp
794 
808  subroutine write_grib2
810  use grib_mod
811 
812  implicit none
813 
814  character(len=1), allocatable :: cgrib(:)
815 
816  integer, parameter :: numcoord = 0
817 
818  integer :: coordlist(numcoord)
819  integer :: lugb, lcgrib, iret
820  integer :: igds(5)
821  integer :: listsec0(2)
822  integer :: listsec1(13)
823  integer :: ideflist, idefnum, ipdsnum, idrsnum
824  integer :: igdstmplen, ipdstmplen, idrstmplen
825  integer :: ipdstmpl(15)
826  integer, allocatable :: igdstmpl(:), idrstmpl(:)
827  integer :: ngrdpts, ibmap, lengrib
828 
829  logical*1, allocatable :: bmap(:), bmap2d(:,:)
830 
831  real, allocatable :: fld(:)
832 
833 !----------------------------------------------------------------------
834 ! Setup variables and arrays required by grib2 library.
835 !----------------------------------------------------------------------
836 
837  call grib2_check(kgds_mdl, igdstmplen)
838 
839  allocate(igdstmpl(igdstmplen))
840 
843  listsec0, listsec1, igds, ipdstmpl, ipdsnum, igdstmpl, &
844  igdstmplen, idefnum, ideflist, ngrdpts)
845 
846  lcgrib = imdl*jmdl*4
847  allocate(cgrib(lcgrib)) ! this variable holds the grib2 message
848 
849  iret=0
850 
851 !----------------------------------------------------------------------
852 ! Create sections 0 and 1. There is no section 2, local use section.
853 !----------------------------------------------------------------------
854 
855  print*,"- CREATE SECTIONS 0 AND 1"
856  call gribcreate(cgrib,lcgrib,listsec0,listsec1,iret)
857  if (iret /= 0) goto 900
858 
859 !----------------------------------------------------------------------
860 ! Create section 3, the grid description section.
861 !----------------------------------------------------------------------
862 
863  print*,"- CREATE SECTION 3"
864  call addgrid(cgrib,lcgrib,igds,igdstmpl,igdstmplen, &
865  ideflist,idefnum,iret)
866  if (iret /= 0) goto 900
867 
868 !----------------------------------------------------------------------
869 ! Setup arrays for section 5, the data representation section.
870 !----------------------------------------------------------------------
871 
872  idrsnum = 0 ! section 5, use simple packing
873  idrstmplen = 5
874  allocate (idrstmpl(idrstmplen))
875  idrstmpl = 0
876  idrstmpl(3) = 0 ! decimal scaling factor
877 
878  allocate(fld(ngrdpts))
879  fld = reshape(snow_cvr_mdl, (/imdl*jmdl/) )
880 
881  ibmap = 0 ! bitmap applies
882  allocate(bmap2d(imdl,jmdl))
883  bmap2d=.true.
884  where (lsmask_mdl < 0.5) bmap2d=.false.
885  allocate(bmap(ngrdpts))
886  bmap = reshape(bmap2d, (/imdl*jmdl/) )
887  deallocate (bmap2d)
888 
889  coordlist=0 ! not used
890 
891 !----------------------------------------------------------------------
892 ! Create section 4 (product definition section) and 5 (data
893 ! representation section) for cover.
894 !----------------------------------------------------------------------
895 
896  ipdstmpl(1) = 1 ! section 4, oct 10; parameter category, moisture
897  ipdstmpl(2) = 42 ! section 4, oct 11; parameter, snow cover
898 
899  print*,"- CREATE SECTIONS 4 AND 5 FOR SNOW COVER"
900  call addfield(cgrib,lcgrib,ipdsnum,ipdstmpl,ipdstmplen, &
901  coordlist,numcoord,idrsnum,idrstmpl, &
902  idrstmplen,fld,ngrdpts,ibmap,bmap,iret)
903  if (iret /= 0) goto 900
904 
905 !----------------------------------------------------------------------
906 ! for regional model, if afwa data not used, don't output a depth
907 ! record. this tells the sfcupdate code to update the first
908 ! guess snow with the new snow cover analysis.
909 !----------------------------------------------------------------------
910 
911  if (kgds_mdl(1) /= 4) then
912  if (.not. use_global_afwa .and. .not. use_nh_afwa .and. .not. use_sh_afwa) goto 88
913  endif
914 
915 !----------------------------------------------------------------------
916 ! Create section 4 (product definition section) and 5 (data
917 ! representation section) for depth.
918 !----------------------------------------------------------------------
919 
920  fld= reshape(snow_dep_mdl, (/imdl*jmdl/) )
921 
922  ipdstmpl(1) = 1 ! section 4, oct 10; parameter category, moisture
923  ipdstmpl(2) = 11 ! section 4, oct 11; parameter, snow depth
924 
925  idrstmpl = 0 ! section 5
926  idrstmpl(3) = 4 ! section 5, decimal scaling factor
927 
928  print*,"- CREATE SECTIONS 4 AND 5 FOR SNOW DEPTH"
929  call addfield(cgrib,lcgrib,ipdsnum,ipdstmpl,ipdstmplen, &
930  coordlist,numcoord,idrsnum,idrstmpl, &
931  idrstmplen,fld,ngrdpts,ibmap,bmap,iret)
932  if (iret /= 0) goto 900
933 
934  88 continue
935 
936 !----------------------------------------------------------------------
937 ! Create section 8 - end section.
938 !----------------------------------------------------------------------
939 
940  call gribend(cgrib,lcgrib,lengrib,iret)
941  if (iret /= 0) goto 900
942 
943 !----------------------------------------------------------------------
944 ! Now output grib message to file.
945 !----------------------------------------------------------------------
946 
947  lugb=53
948  print*,"- OPEN OUTPUT GRIB FILE ", trim(model_snow_file)
949  call baopenw(lugb, model_snow_file, iret)
950 
951  if (iret /= 0) then
952  print*,'- FATAL ERROR: BAD OPEN OF OUTPUT GRIB FILE. IRET IS ', iret
953  call w3tage('SNOW2MDL')
954  call errexit(49)
955  end if
956 
957  print*,'- WRITE OUTPUT GRIB FILE.'
958  call wryte(lugb, lengrib, cgrib)
959 
960  call baclose (lugb, iret)
961 
962  deallocate(fld, bmap, idrstmpl, igdstmpl, cgrib)
963 
964  return
965 
966  900 continue
967  print*,'- FATAL ERROR CREATING GRIB2 MESSAGE. IRET IS ', iret
968  call w3tage('SNOW2MDL')
969  call errexit(48)
970 
971  end subroutine write_grib2
972 
988  subroutine write_grib1
990  implicit none
991 
992  integer :: iret
993  integer, parameter :: lugb = 64 ! unit number of output grib file
994  integer :: kpds(200)
995 
996  logical*1 :: lbms(imdl,jmdl)
997 
998 !----------------------------------------------------------------------
999 ! set up pds section. don't need to set the gds section.
1000 ! since the model grid is not changing, use the kgds array
1001 ! already determined in module model_grid.
1002 !----------------------------------------------------------------------
1003 
1004  kpds = 0
1005 
1006  kpds(1) = 7 ! center id
1007  kpds(2) = 25 ! process id number. this determined from the
1008  ! input data as we are simply interpolating
1009  ! that data to a different grid. should
1010  ! i request a process id for my codes?
1011  kpds(3) = grid_id_mdl ! grid specified in gds
1012  kpds(4) = 192 ! include gds and a bit map section
1013  kpds(5) = 238 ! parameter number for snow cover
1014  kpds(6) = 1 ! level - ground or water surface
1015  kpds(7) = 0 ! height pressure of level
1016  kpds(8) = grib_year ! year of century the time info is determined
1017  kpds(9) = grib_month ! month by operational requirements
1018  kpds(10) = grib_day ! day
1019  kpds(11) = grib_hour ! hour
1020  kpds(12) = 0 ! minute
1021  kpds(13) = 1 ! fcst time unit - hour
1022  kpds(14) = 0 ! period of time, p1. set to '0' for analysis
1023  kpds(15) = 0 ! number of time units, p2.
1024  kpds(16) = 1 ! initialized analysis product
1025  kpds(17) = 0 ! number in average
1026  kpds(18) = 1 ! grib edition 1
1027  kpds(19) = 3 ! parameter table version number
1028  kpds(20) = 0 ! number missing from avg/accum
1029  kpds(21) = grib_century ! century - set as in the input file
1030  kpds(22) = 0 ! decimal scale factor
1031  kpds(23) = 4 ! subcenter - emc
1032  kpds(24) = 0 ! reserved
1033  kpds(25) = 0 ! reserved
1034 
1035  lbms = .false. ! set bitmap section
1036 
1037  where(lsmask_mdl > 0.0) lbms = .true.
1038 
1039  print*,"- OPEN OUTPUT GRIB FILE ", trim(model_snow_file)
1040  call baopenw(lugb, model_snow_file, iret)
1041 
1042  if (iret /= 0) then
1043  print*,'- FATAL ERROR OPENING OUTPUT GRIB FILE. IRET IS ', iret
1044  call w3tage('SNOW2MDL')
1045  call errexit(59)
1046  end if
1047 
1048  print*,"- WRITE OUTPUT GRIB FILE ", trim(model_snow_file)
1049  call putgb (lugb, (imdl*jmdl), kpds, kgds_mdl, lbms, &
1050  snow_cvr_mdl, iret)
1051 
1052  if (iret /= 0) then
1053  print*,'- FATAL ERROR WRITING OUTPUT GRIB FILE. IRET IS ', iret
1054  call w3tage('SNOW2MDL')
1055  call errexit(58)
1056  end if
1057 
1058 ! for regional model, if afwa data not used, don't output a depth
1059 ! record. this tells the sfcupdate code to update the first
1060 ! guess snow with the new snow cover analysis.
1061 
1062  if (kgds_mdl(1) /= 4) then
1063  if (.not. use_global_afwa .and. .not. use_nh_afwa .and. .not. use_sh_afwa) goto 88
1064  endif
1065 
1066  kpds(5) = 66 ! parameter number for snow depth
1067  kpds(22) = 4 ! scaling factor. to nearest mm
1068 
1069  call putgb (lugb, (imdl*jmdl), kpds, kgds_mdl, lbms, &
1070  snow_dep_mdl, iret)
1071 
1072  if (iret /= 0) then
1073  print*,'- FATAL ERROR WRITING OUTPUT GRIB FILE. IRET IS ', iret
1074  call w3tage('SNOW2MDL')
1075  call errexit(57)
1076  end if
1077 
1078  88 call baclose(lugb, iret)
1079 
1080  return
1081 
1082  end subroutine write_grib1
1083 
1100  subroutine uninterpred(iord,kmsk,fi,f,lonl,latd,len,lonsperlat)
1102  implicit none
1103 
1104  integer, intent(in) :: len
1105  integer, intent(in) :: iord
1106  integer, intent(in) :: lonl
1107  integer, intent(in) :: latd
1108  integer, intent(in) :: lonsperlat(latd/2)
1109  integer, intent(in) :: kmsk(lonl*latd)
1110  integer :: j,lons,jj,latd2,ii,i
1111 
1112  real, intent(in) :: fi(len)
1113  real, intent(out) :: f(lonl,latd)
1114 
1115  latd2 = latd / 2
1116  ii = 1
1117 
1118  do j=1,latd
1119 
1120  jj = j
1121  if (j .gt. latd2) jj = latd - j + 1
1122  lons=lonsperlat(jj)
1123 
1124  if(lons.ne.lonl) then
1125  call intlon(iord,1,1,lons,lonl,kmsk(ii),fi(ii),f(1,j))
1126  else
1127  do i=1,lonl
1128  f(i,j) = fi(ii+i-1)
1129  enddo
1130  endif
1131 
1132  ii = ii + lons
1133 
1134  enddo
1135 
1136  end subroutine uninterpred
1137 
1151  subroutine intlon(iord,imon,imsk,m1,m2,k1,f1,f2)
1153  implicit none
1154 
1155  integer,intent(in) :: iord,imon,imsk,m1,m2
1156  integer,intent(in) :: k1(m1)
1157  integer :: i2,in,il,ir
1158 
1159  real,intent(in) :: f1(m1)
1160  real,intent(out) :: f2(m2)
1161  real :: r,x1
1162 
1163  r=real(m1)/real(m2)
1164  do i2=1,m2
1165  x1=(i2-1)*r
1166  il=int(x1)+1
1167  ir=mod(il,m1)+1
1168  if(iord.eq.2.and.(imsk.eq.0.or.k1(il).eq.k1(ir))) then
1169  f2(i2)=f1(il)*(il-x1)+f1(ir)*(x1-il+1)
1170  else
1171  in=mod(nint(x1),m1)+1
1172  f2(i2)=f1(in)
1173  endif
1174  enddo
1175 
1176  end subroutine intlon
1177 
1178  end module snow2mdl
real, dimension(:,:), allocatable snow_dep_afwa_sh
Southern hemisphere afwa snow depth.
Definition: snowdat.F90:90
This module reads in data from the program&#39;s configuration namelist.
logical bad_afwa_sh
When true, the southern hemisphere afwa data failed its quality control check.
Definition: snowdat.F90:62
real, dimension(:,:), allocatable snow_cvr_mdl
snow cover on model grid in percent
Definition: snow2mdl.F90:85
integer, dimension(200) kgds_afwa_nh
grib1 grid description section for northern hemisphere 16th mesh afwa data.
Definition: snowdat.F90:48
integer mesh_nesdis
nesdis/ims data is 96th mesh (or bediant)
Definition: snowdat.F90:58
real, public snow_cvr_threshold
if percent coverage according to nesdis/ims or autosnow exceeds this value, then non-zero snow depth ...
real, public min_snow_depth
minimum snow depth in meters at model points with coverage exceeding threshold.
real lonlast
Corner point longitude (imdl,jmdl) of model grid.
Definition: model_grid.F90:45
subroutine init_grib2(century, year, month, day, hour, kgds, lat11, latlast, lon11, lonlast, listsec0, listsec1, igds, ipdstmpl, ipdsnum, igdstmpl, igdstmplen, idefnum, ideflist, ngrdpts)
Initialize grib2 arrays required by the ncep g2 library according to grib1 gds information.
Definition: grib_utils.F90:401
subroutine grib2_check(kgds, igdstmplen)
Determine length of grib2 gds template array, which is a function of the map projection.
Definition: grib_utils.F90:351
real, dimension(:,:), allocatable snow_cvr_autosnow
autosnow snow cover flag (0-no, 100-yes)
Definition: snowdat.F90:87
integer *1, dimension(:,:), allocatable sea_ice_nesdis
nesdis/ims sea ice flag (0-open water, 1-ice)
Definition: snowdat.F90:59
real, dimension(:,:), allocatable snow_dep_afwa_global
The global afwa snow depth.
Definition: snowdat.F90:88
integer, public grib_century
date of the final merged snow product that will be placed in grib header.
integer jmdl
j-dimension of model grid
Definition: model_grid.F90:28
integer, dimension(200) kgds_nesdis
nesdis/ims grid description section (grib section 2)
Definition: snowdat.F90:57
logical *1, dimension(:,:), allocatable bitmap_nesdis
nesdis data grib bitmap (false-non land, true-land).
Definition: snowdat.F90:74
integer, dimension(:), allocatable lonsperlat_mdl
Number of longitudes (i-points) for each latitude (row).
Definition: model_grid.F90:34
logical thinned
When true, global grids will run thinned (number of i points decrease toward pole) ...
Definition: model_grid.F90:38
real afwa_res
Resolution of afwa data in km.
Definition: snowdat.F90:84
real, dimension(:), allocatable lats_mdl
Latitudes of model grid points.
Definition: model_grid.F90:41
integer, dimension(200) kgds_afwa_global
grib1 grid description section for global afwa data.
Definition: snowdat.F90:46
Read in data defining the model grid.
Definition: model_grid.F90:19
integer, dimension(200) kgds_autosnow
autosnow grid description section (grib section 2)
Definition: snowdat.F90:56
integer jafwa
j-dimension of afwa grid
Definition: snowdat.F90:43
integer iautosnow
i-dimension of autosnow grid
Definition: snowdat.F90:41
Read and qc afwa, nesdis/ims and autosnow snow data.
Definition: snowdat.F90:26
integer iafwa
i-dimension of afwa grid
Definition: snowdat.F90:40
real, dimension(:,:), allocatable lsmask_mdl
land mask of model grid (0 - non land, 1-land) for global grids run thinned, will contain a modified ...
Definition: model_grid.F90:47
subroutine intlon(iord, imon, imsk, m1, m2, k1, f1, f2)
Convert data from the thinned (or reduced) to the full grid along a single row.
Definition: snow2mdl.F90:1152
subroutine uninterpred(iord, kmsk, fi, f, lonl, latd, len, lonsperlat)
Fills out full grid using thinned grid data.
Definition: snow2mdl.F90:1101
logical *1, dimension(:,:), allocatable bitmap_autosnow
autosnow data grib bitmap (false-non land, true-land).
Definition: snowdat.F90:75
real nesdis_res
Resolution of the nesdis data in km.
Definition: snowdat.F90:85
integer, public grib_month
date of the final merged snow product that will be placed in grib header.
logical *1, dimension(:,:), allocatable bitmap_afwa_sh
The southern hemisphere afwa data grib bitmap.
Definition: snowdat.F90:72
real, dimension(:), allocatable lons_mdl
longitudes of model grid points
Definition: model_grid.F90:46
real, dimension(:,:), allocatable snow_dep_afwa_nh
Northern hemisphere afwa snow depth.
Definition: snowdat.F90:89
real resol_mdl
approximate model resolution in km.
Definition: model_grid.F90:52
logical use_sh_afwa
True if southern hemisphere afwa data to be used.
Definition: snowdat.F90:78
Interpolate snow data to model grid and grib the result.
Definition: snow2mdl.F90:20
subroutine write_grib2
Write grib2 snow cover and depth on the model grid.
Definition: snow2mdl.F90:809
integer inesdis
i-dimension of nesdis grid
Definition: snowdat.F90:42
integer imdl
i-dimension of model grid
Definition: model_grid.F90:27
integer, public grib_day
date of the final merged snow product that will be placed in grib header.
integer jautosnow
j-dimension of autosnow grid
Definition: snowdat.F90:44
real, public lat_threshold
equatorward of this latitude, model points with undefined cover or depth (because the interpolation r...
logical use_nh_afwa
True if northern hemisphere afwa data to be used.
Definition: snowdat.F90:77
logical use_nesdis
True if nesdis/ims data to be used.
Definition: snowdat.F90:81
logical, public output_grib2
when true, output model snow analysis is grib 2.
subroutine write_grib1
Write grib1 snow cover and depth on the model grid.
Definition: snow2mdl.F90:989
real lon11
Corner point longitude (1,1) of model grid.
Definition: model_grid.F90:44
integer ijmdl
total number of model land points
Definition: model_grid.F90:29
logical use_autosnow
True if autosnow data to be used.
Definition: snowdat.F90:80
logical *1, dimension(:,:), allocatable bitmap_afwa_global
The global afwa data grib bitmap.
Definition: snowdat.F90:68
character *200, public model_snow_file
path/name nesdis/ims snow cover
integer, dimension(:), allocatable jpts_mdl
j index of point on full grid
Definition: model_grid.F90:31
integer, dimension(200) kgds_mdl
holds grib gds info of model grid
Definition: model_grid.F90:33
integer, public grib_year
date of the final merged snow product that will be placed in grib header.
integer jnesdis
j-dimension of nesdis grid
Definition: snowdat.F90:45
real latlast
Corner point latitude (imdl,jmdl) of model grid.
Definition: model_grid.F90:43
integer, dimension(:), allocatable ipts_mdl
i index of point on full grid
Definition: model_grid.F90:30
integer, public grib_hour
date of the final merged snow product that will be placed in grib header.
logical bad_afwa_nh
When true, the northern hemisphere afwa data failed its quality control check.
Definition: snowdat.F90:60
integer, dimension(200) kgds_afwa_sh
grib1 grid description section for southern hemisphere 16th mesh afwa data.
Definition: snowdat.F90:52
real autosnow_res
Resolution of autosnow in km.
Definition: snowdat.F90:83
logical *1, dimension(:,:), allocatable bitmap_afwa_nh
The northern hemisphere afwa data grib bitmap.
Definition: snowdat.F90:70
real, dimension(:,:), allocatable snow_dep_mdl
snow depth on model grid in meters
Definition: snow2mdl.F90:86
integer grid_id_mdl
grib id of model grid, 4-gaussian, 203-egrid
Definition: model_grid.F90:26
logical use_global_afwa
True if global hemisphere afwa data to be used.
Definition: snowdat.F90:79
real lat11
Corner point latitude (1,1) of model grid.
Definition: model_grid.F90:42
real, dimension(:,:), allocatable lsmask_mdl_sav
saved copy of land mask of model grid (0 - non land, 1-land) only used for global thinned grids...
Definition: model_grid.F90:50
real, dimension(:,:), allocatable snow_cvr_nesdis
nesdis/ims snow cover flag (0-no, 100-yes)
Definition: snowdat.F90:86
subroutine, public interp
Interpolate snow data to model grid.
Definition: snow2mdl.F90:158