GFS Operational Physics Documentation  gsm/branches/DTC/phys-doc-all phys-doc-all R82971
radiation_surface.f
Go to the documentation of this file.
1 
4 
5 ! ========================================================== !!!!!
6 ! 'module_radiation_surface' description !!!!!
7 ! ========================================================== !!!!!
8 ! !
9 ! this module sets up surface albedo for sw radiation and surface !
10 ! emissivity for lw radiation. !
11 ! !
12 ! !
13 ! in the module, the externally callabe subroutines are : !
14 ! !
15 ! 'sfc_init' -- initialization radiation surface data !
16 ! inputs: !
17 ! ( me ) !
18 ! outputs: !
19 ! (none) !
20 ! !
21 ! 'setalb' -- set up four-component surface albedoes !
22 ! inputs: !
23 ! (slmsk,snowf,sncovr,snoalb,zorlf,coszf,tsknf,tairf,hprif, !
24 ! alvsf,alnsf,alvwf,alnwf,facsf,facwf,fice,tisfc !
25 ! IMAX) !
26 ! outputs: !
27 ! (sfcalb) !
28 ! !
29 ! 'setemis' -- set up surface emissivity for lw radiation !
30 ! inputs: !
31 ! (xlon,xlat,slmsk,snowf,sncovr,zorlf,tsknf,tairf,hprif, !
32 ! IMAX) !
33 ! outputs: !
34 ! (sfcemis) !
35 ! !
36 ! external modules referenced: !
37 ! !
38 ! 'module machine' in 'machine.f' !
39 ! 'module physcons' in 'physcons.f' !
40 ! 'module module_iounitdef' in 'iounitdef.f' !
41 ! !
42 ! !
43 ! program history log: !
44 ! 1995 y.t. hou - created albaer.f (include albedo !
45 ! and aerosols calculations) !
46 ! nov 1997 y.t. hou - modified snow albedo !
47 ! jan 1998 y.t. hou - included grumbine's sea-ice scheme !
48 ! feb 1998 h.l. pan - seasonal interpolation in cycle !
49 ! mar 2000 y.t. hou - modified to use opac aerosol data !
50 ! apr 2003 y.t. hou - seperate albedo and aerosols into !
51 ! two subroutines, rewritten in f90 modulized form !
52 ! jan 2005 s. moorthi - xingren's sea-ice fraction added !
53 ! apr 2005 y.t. hou - revised module structure !
54 ! feb 2006 y.t. hou - add varying surface emissivity, !
55 ! modified sfc albedo structure for modis shceme !
56 ! Mar 2006 s. moorthi - added surface temp over ice fraction !
57 ! mar 2007 c. marshall & h. wei !
58 ! - added modis based sfc albedo scheme !
59 ! may 2007 y. hou & s. moorthi !
60 ! - fix bug in modis albedo over ocean !
61 ! aug 2007 h. wei & s. moorthi !
62 ! - fix bug in modis albedo over sea-ice !
63 ! aug 2007 y. hou - fix bug in emissivity over ocean in !
64 ! the modis scheme option !
65 ! dec 2008 f. yang - modified zenith angle dependence on !
66 ! surface albedo over land. (2008 jamc)!
67 ! aug 2012 y. hou - minor modification in initialization !
68 ! subr 'sfc_init'. !
69 ! nov 2012 y. hou - modified control parameters through !
70 ! module 'physparam'. !
71 ! !
72 !!!!! ========================================================== !!!!!
73 !!!!! end descriptions !!!!!
74 !!!!! ========================================================== !!!!!
75 
76 
83 !========================================!
84  module module_radiation_surface !
85 !........................................!
86 !
87  use physparam, only : ialbflg, iemsflg, semis_file, &
88  & kind_phys
89  use physcons, only : con_t0c, con_ttp, con_pi, con_tice
90  use module_iounitdef, only : niradsf
91 !
92  implicit none
93 !
94  private
95 
96 ! --- version tag and last revision date
97  character(40), parameter :: &
98  & VTAGSFC='NCEP-Radiation_surface v5.1 Nov 2012 '
99 ! & VTAGSFC='NCEP-Radiation_surface v5.0 Aug 2012 '
100 
101 ! --- constant parameters
103  integer, parameter, public :: nf_albd = 4
104 
106  integer, parameter, public :: imxems = 360
107 
109  integer, parameter, public :: jmxems = 180
110 
111  real (kind=kind_phys), parameter :: f_zero = 0.0
112  real (kind=kind_phys), parameter :: f_one = 1.0
113  real (kind=kind_phys), parameter :: rad2dg= 180.0 / con_pi
114 
116  integer, allocatable :: idxems(:,:)
118  integer :: iemslw = 0
119 !
120  public sfc_init, setalb, setemis
121 
122 ! =================
123  contains
124 ! =================
125 
126 
132 !-----------------------------------
133  subroutine sfc_init &
134  & ( me )! --- inputs:
135 ! --- outputs: ( none )
136 
137 ! =================================================================== !
138 ! !
139 ! this program is the initialization program for surface radiation !
140 ! related quantities (albedo, emissivity, etc.) !
141 ! !
142 ! usage: call sfc_init !
143 ! !
144 ! subprograms called: none !
145 ! !
146 ! ==================== defination of variables ==================== !
147 ! !
148 ! inputs: !
149 ! me - print control flag !
150 ! !
151 ! outputs: (none) to module variables only !
152 ! !
153 ! external module variables: !
154 ! ialbflg - control flag for surface albedo schemes !
155 ! =0: climatology, based on surface veg types !
156 ! =1: !
157 ! iemsflg - control flag for sfc emissivity schemes (ab:2-dig)!
158 ! a:=0 set sfc air/ground t same for lw radiation !
159 ! =1 set sfc air/ground t diff for lw radiation !
160 ! b:=0 use fixed sfc emissivity=1.0 (black-body) !
161 ! =1 use varying climtology sfc emiss (veg based) !
162 ! !
163 ! ==================== end of description ===================== !
164 !
165  implicit none
166 
167 ! --- inputs:
168  integer, intent(in) :: me
169 
170 ! --- outputs: ( none )
171 
172 ! --- locals:
173  integer :: i, k
174 ! integer :: ia, ja
175  logical :: file_exist
176  character :: cline*80
177 !
178 !===> ... begin here
179 !
180  if ( me == 0 ) print *, vtagsfc ! print out version tag
181 
186 
187  if ( ialbflg == 0 ) then
188 
189  if ( me == 0 ) then
190  print *,' - Using climatology surface albedo scheme for sw'
191  endif
192 
193  else if ( ialbflg == 1 ) then
194 
195  if ( me == 0 ) then
196  print *,' - Using MODIS based land surface albedo for sw'
197  endif
198 
199  else
200  print *,' !! ERROR in Albedo Scheme Setting, IALB=',ialbflg
201  stop
202  endif ! end if_ialbflg_block
203 
208 
209  iemslw = mod(iemsflg, 10) ! emissivity control
210  if ( iemslw == 0 ) then ! fixed sfc emis at 1.0
211 
212  if ( me == 0 ) then
213  print *,' - Using Fixed Surface Emissivity = 1.0 for lw'
214  endif
215 
216  elseif ( iemslw == 1 ) then ! input sfc emiss type map
217 
218 ! --- allocate data space
219  if ( .not. allocated(idxems) ) then
220  allocate ( idxems(imxems,jmxems) )
221  endif
222 
223 ! --- check to see if requested emissivity data file existed
224 
225  inquire (file=semis_file, exist=file_exist)
226 
227  if ( .not. file_exist ) then
228  if ( me == 0 ) then
229  print *,' - Using Varying Surface Emissivity for lw'
230  print *,' Requested data file "',semis_file,'" not found!'
231  print *,' Change to fixed surface emissivity = 1.0 !'
232  endif
233 
234  iemslw = 0
235  else
236  close(niradsf)
237  open (niradsf,file=semis_file,form='formatted',status='old')
238  rewind niradsf
239 
240  read (niradsf,12) cline
241  12 format(a80)
242 
243  read (niradsf,14) idxems
244  14 format(80i1)
245 
246  if ( me == 0 ) then
247  print *,' - Using Varying Surface Emissivity for lw'
248  print *,' Opened data file: ',semis_file
249  print *, cline
250 !check print *,' CHECK: Sample emissivity index data'
251 ! ia = IMXEMS / 5
252 ! ja = JMXEMS / 5
253 ! print *, idxems(1:IMXEMS:ia,1:JMXEMS:ja)
254  endif
255 
256  close(niradsf)
257  endif ! end if_file_exist_block
258 
259  else
260  print *,' !! ERROR in Emissivity Scheme Setting, IEMS=',iemsflg
261  stop
262  endif ! end if_iemslw_block
263 
264 !
265  return
266 !...................................
267  end subroutine sfc_init
268 !-----------------------------------
269 !! @}
270 
271 
307 !-----------------------------------
308  subroutine setalb &
309  & ( slmsk,snowf,sncovr,snoalb,zorlf,coszf,tsknf,tairf,hprif, & ! --- inputs:
310  & alvsf,alnsf,alvwf,alnwf,facsf,facwf,fice,tisfc, &
311  & imax, &
312  & sfcalb & ! --- outputs:
313  & )
315 ! =================================================================== !
316 ! !
317 ! this program computes four components of surface albedos (i.e. !
318 ! vis-nir, direct-diffused) according to controflag ialbflg. !
319 ! 1) climatological surface albedo scheme (briegleb 1992) !
320 ! 2) modis retrieval based scheme from boston univ. !
321 ! !
322 ! !
323 ! usage: call setalb !
324 ! !
325 ! subprograms called: none !
326 ! !
327 ! ==================== defination of variables ==================== !
328 ! !
329 ! inputs: !
330 ! slmsk (IMAX) - sea(0),land(1),ice(2) mask on fcst model grid !
331 ! snowf (IMAX) - snow depth water equivalent in mm !
332 ! sncovr(IMAX) - ialgflg=0: not used !
333 ! ialgflg=1: snow cover over land in fraction !
334 ! snoalb(IMAX) - ialbflg=0: not used !
335 ! ialgflg=1: max snow albedo over land in fraction !
336 ! zorlf (IMAX) - surface roughness in cm !
337 ! coszf (IMAX) - cosin of solar zenith angle !
338 ! tsknf (IMAX) - ground surface temperature in k !
339 ! tairf (IMAX) - lowest model layer air temperature in k !
340 ! hprif (IMAX) - topographic sdv in m !
341 ! --- for ialbflg=0 climtological albedo scheme --- !
342 ! alvsf (IMAX) - 60 degree vis albedo with strong cosz dependency !
343 ! alnsf (IMAX) - 60 degree nir albedo with strong cosz dependency !
344 ! alvwf (IMAX) - 60 degree vis albedo with weak cosz dependency !
345 ! alnwf (IMAX) - 60 degree nir albedo with weak cosz dependency !
346 ! --- for ialbflg=1 modis based land albedo scheme --- !
347 ! alvsf (IMAX) - visible black sky albedo at zenith 60 degree !
348 ! alnsf (IMAX) - near-ir black sky albedo at zenith 60 degree !
349 ! alvwf (IMAX) - visible white sky albedo !
350 ! alnwf (IMAX) - near-ir white sky albedo !
351 ! !
352 ! facsf (IMAX) - fractional coverage with strong cosz dependency !
353 ! facwf (IMAX) - fractional coverage with weak cosz dependency !
354 ! fice (IMAX) - sea-ice fraction !
355 ! tisfc (IMAX) - sea-ice surface temperature !
356 ! IMAX - array horizontal dimension !
357 ! !
358 ! outputs: !
359 ! sfcalb(IMAX,NF_ALBD) !
360 ! ( :, 1) - near ir direct beam albedo !
361 ! ( :, 2) - near ir diffused albedo !
362 ! ( :, 3) - uv+vis direct beam albedo !
363 ! ( :, 4) - uv+vis diffused albedo !
364 ! !
365 ! module internal control variables: !
366 ! ialbflg - =0 use the default climatology surface albedo !
367 ! =1 use modis retrieved albedo and input snow cover!
368 ! for land areas !
369 ! !
370 ! ==================== end of description ===================== !
371 !
372  implicit none
373 
374 ! --- inputs
375  integer, intent(in) :: IMAX
376 
377  real (kind=kind_phys), dimension(:), intent(in) :: &
378  & slmsk, snowf, zorlf, coszf, tsknf, tairf, hprif, &
379  & alvsf, alnsf, alvwf, alnwf, facsf, facwf, fice, tisfc, &
380  & sncovr, snoalb
381 
382 ! --- outputs
383  real (kind=kind_phys), dimension(IMAX,NF_ALBD), intent(out) :: &
384  & sfcalb
385 ! real (kind=kind_phys), dimension(:,:), intent(out) :: sfcalb
386 
387 ! --- locals:
388  real (kind=kind_phys) :: asnvb, asnnb, asnvd, asnnd, asevb &
389  &, asenb, asevd, asend, fsno, fsea, rfcs, rfcw, flnd &
390  &, asnow, argh, hrgh, fsno0, fsno1, flnd0, fsea0, csnow &
391  &, a1, a2, b1, b2, b3, ab1bm, ab2bm
392 
393  real (kind=kind_phys) ffw, dtgd
394 
395  integer :: i, k
396 
397 !
398 !===> ... begin here
399 !
400 
402  if ( ialbflg == 0 ) then ! use climatological albedo scheme
403 
404  do i = 1, imax
405 
408 
409  asnow = 0.02*snowf(i)
410  argh = min(0.50, max(.025, 0.01*zorlf(i)))
411  hrgh = min(f_one, max(0.20, 1.0577-1.1538e-3*hprif(i) ) )
412  fsno0 = asnow / (argh + asnow) * hrgh
413  if (nint(slmsk(i))==0 .and. tsknf(i)>con_tice) fsno0 = f_zero
414  fsno1 = f_one - fsno0
415  flnd0 = min(f_one, facsf(i)+facwf(i))
416  fsea0 = max(f_zero, f_one-flnd0)
417  fsno = fsno0
418  fsea = fsea0 * fsno1
419  flnd = flnd0 * fsno1
420 
422 
423  if (tsknf(i) >= 271.5) then
424  asevd = 0.06
425  asend = 0.06
426  elseif (tsknf(i) < 271.1) then
427  asevd = 0.70
428  asend = 0.65
429  else
430  a1 = (tsknf(i) - 271.1)**2
431  asevd = 0.7 - 4.0*a1
432  asend = 0.65 - 3.6875*a1
433  endif
434 
436 
437  if (nint(slmsk(i)) == 2) then
438  ffw = f_one - fice(i)
439  if (ffw < f_one) then
440  dtgd = max(f_zero, min(5.0, (con_ttp-tisfc(i)) ))
441  b1 = 0.03 * dtgd
442  else
443  b1 = f_zero
444  endif
445 
446  b3 = 0.06 * ffw
447  asnvd = (0.70 + b1) * fice(i) + b3
448  asnnd = (0.60 + b1) * fice(i) + b3
449  asevd = 0.70 * fice(i) + b3
450  asend = 0.60 * fice(i) + b3
451  else
452  asnvd = 0.90
453  asnnd = 0.75
454  endif
455 
457 
458  if (coszf(i) < 0.5) then
459  csnow = 0.5 * (3.0 / (f_one+4.0*coszf(i)) - f_one)
460  asnvb = min( 0.98, asnvd+(1.0-asnvd)*csnow )
461  asnnb = min( 0.98, asnnd+(1.0-asnnd)*csnow )
462  else
463  asnvb = asnvd
464  asnnb = asnnd
465  endif
466 
468 
469  if (coszf(i) > 0.0001) then
470 ! rfcs = 1.4 / (f_one + 0.8*coszf(i))
471 ! rfcw = 1.3 / (f_one + 0.6*coszf(i))
472  rfcs = 2.14 / (f_one + 1.48*coszf(i))
473  rfcw = rfcs
474 
475  if (tsknf(i) >= con_t0c) then
476  asevb = max(asevd, 0.026/(coszf(i)**1.7+0.065) &
477  & + 0.15 * (coszf(i)-0.1) * (coszf(i)-0.5) &
478  & * (coszf(i)-f_one))
479  asenb = asevb
480  else
481  asevb = asevd
482  asenb = asend
483  endif
484  else
485  rfcs = f_one
486  rfcw = f_one
487  asevb = asevd
488  asenb = asend
489  endif
490 
491  a1 = alvsf(i) * facsf(i)
492  b1 = alvwf(i) * facwf(i)
493  a2 = alnsf(i) * facsf(i)
494  b2 = alnwf(i) * facwf(i)
495  ab1bm = a1*rfcs + b1*rfcw
496  ab2bm = a2*rfcs + b2*rfcw
497  sfcalb(i,1) = min(0.99, ab2bm) *flnd + asenb*fsea + asnnb*fsno
498  sfcalb(i,2) = (a2 + b2) * 0.96 *flnd + asend*fsea + asnnd*fsno
499  sfcalb(i,3) = min(0.99, ab1bm) *flnd + asevb*fsea + asnvb*fsno
500  sfcalb(i,4) = (a1 + b1) * 0.96 *flnd + asevd*fsea + asnvd*fsno
501 
502  enddo ! end_do_i_loop
503 
505  else
506 
507  do i = 1, imax
508 
511 
512  fsno0 = sncovr(i)
513 
514  if (nint(slmsk(i))==0 .and. tsknf(i)>con_tice) fsno0 = f_zero
515 
516  if (nint(slmsk(i)) == 2) then
517  asnow = 0.02*snowf(i)
518  argh = min(0.50, max(.025, 0.01*zorlf(i)))
519  hrgh = min(f_one, max(0.20, 1.0577-1.1538e-3*hprif(i) ) )
520  fsno0 = asnow / (argh + asnow) * hrgh
521  endif
522 
523  fsno1 = f_one - fsno0
524  flnd0 = min(f_one, facsf(i)+facwf(i))
525  fsea0 = max(f_zero, f_one-flnd0)
526  fsno = fsno0
527  fsea = fsea0 * fsno1
528  flnd = flnd0 * fsno1
529 
531 
532  if (tsknf(i) >= 271.5) then
533  asevd = 0.06
534  asend = 0.06
535  elseif (tsknf(i) < 271.1) then
536  asevd = 0.70
537  asend = 0.65
538  else
539  a1 = (tsknf(i) - 271.1)**2
540  asevd = 0.7 - 4.0*a1
541  asend = 0.65 - 3.6875*a1
542  endif
543 
546 
547  if (nint(slmsk(i)) == 2) then
548  ffw = f_one - fice(i)
549  if (ffw < f_one) then
550  dtgd = max(f_zero, min(5.0, (con_ttp-tisfc(i)) ))
551  b1 = 0.03 * dtgd
552  else
553  b1 = f_zero
554  endif
555 
556  b3 = 0.06 * ffw
557  asnvd = (0.70 + b1) * fice(i) + b3
558  asnnd = (0.60 + b1) * fice(i) + b3
559  asevd = 0.70 * fice(i) + b3
560  asend = 0.60 * fice(i) + b3
561  else
562  asnvd = snoalb(i)
563  asnnd = snoalb(i)
564  endif
565 
567 
568  if (nint(slmsk(i)) == 2) then
569  if (coszf(i) < 0.5) then
570  csnow = 0.5 * (3.0 / (f_one+4.0*coszf(i)) - f_one)
571  asnvb = min( 0.98, asnvd+(f_one-asnvd)*csnow )
572  asnnb = min( 0.98, asnnd+(f_one-asnnd)*csnow )
573  else
574  asnvb = asnvd
575  asnnb = asnnd
576  endif
577  else
578  asnvb = snoalb(i)
579  asnnb = snoalb(i)
580  endif
581 
584 
585  if (coszf(i) > 0.0001) then
586 
587 ! rfcs = 1.89 - 3.34*coszf(i) + 4.13*coszf(i)*coszf(i) &
588 ! & - 2.02*coszf(i)*coszf(i)*coszf(i)
589  rfcs = 1.775/(1.0+1.55*coszf(i))
590 
591  if (tsknf(i) >= con_t0c) then
592  asevb = max(asevd, 0.026/(coszf(i)**1.7+0.065) &
593  & + 0.15 * (coszf(i)-0.1) * (coszf(i)-0.5) &
594  & * (coszf(i)-f_one))
595  asenb = asevb
596  else
597  asevb = asevd
598  asenb = asend
599  endif
600  else
601  rfcs = f_one
602  asevb = asevd
603  asenb = asend
604  endif
605 
606  ab1bm = min(0.99, alnsf(i)*rfcs)
607  ab2bm = min(0.99, alvsf(i)*rfcs)
608  sfcalb(i,1) = ab1bm *flnd + asenb*fsea + asnnb*fsno
609  sfcalb(i,2) = alnwf(i) *flnd + asend*fsea + asnnd*fsno
610  sfcalb(i,3) = ab2bm *flnd + asevb*fsea + asnvb*fsno
611  sfcalb(i,4) = alvwf(i) *flnd + asevd*fsea + asnvd*fsno
612 
613  enddo ! end_do_i_loop
614 
615  endif ! end if_ialbflg
616 !
617  return
618 !...................................
619  end subroutine setalb
620 !-----------------------------------
621 !! @}
622 
639 !-----------------------------------
640  subroutine setemis &
641  & ( xlon,xlat,slmsk,snowf,sncovr,zorlf,tsknf,tairf,hprif, & ! --- inputs:
642  & imax, &
643  & sfcemis & ! --- outputs:
644  & )
646 ! =================================================================== !
647 ! !
648 ! this program computes surface emissivity for lw radiation. !
649 ! !
650 ! usage: call setemis !
651 ! !
652 ! subprograms called: none !
653 ! !
654 ! ==================== defination of variables ==================== !
655 ! !
656 ! inputs: !
657 ! xlon (IMAX) - longitude in radiance, ok for both 0->2pi or !
658 ! -pi -> +pi ranges !
659 ! xlat (IMAX) - latitude in radiance, default to pi/2 -> -pi/2 !
660 ! range, otherwise see in-line comment !
661 ! slmsk (IMAX) - sea(0),land(1),ice(2) mask on fcst model grid !
662 ! snowf (IMAX) - snow depth water equivalent in mm !
663 ! sncovr(IMAX) - ialbflg=1: snow cover over land in fraction !
664 ! zorlf (IMAX) - surface roughness in cm !
665 ! tsknf (IMAX) - ground surface temperature in k !
666 ! tairf (IMAX) - lowest model layer air temperature in k !
667 ! hprif (IMAX) - topographic sdv in m !
668 ! IMAX - array horizontal dimension !
669 ! !
670 ! outputs: !
671 ! sfcemis(IMAX) - surface emissivity !
672 ! !
673 ! ------------------------------------------------------------------- !
674 ! !
675 ! surface type definations: !
676 ! 1. open water 2. grass/wood/shrub land !
677 ! 3. tundra/bare soil 4. sandy desert !
678 ! 5. rocky desert 6. forest !
679 ! 7. ice 8. snow !
680 ! !
681 ! input index data lon from 0 towards east, lat from n to s !
682 ! !
683 ! ==================== end of description ===================== !
684 !
685  implicit none
686 
687 ! --- inputs
688  integer, intent(in) :: IMAX
689 
690  real (kind=kind_phys), dimension(:), intent(in) :: &
691  & xlon,xlat, slmsk, snowf,sncovr, zorlf, tsknf, tairf, hprif
692 
693 ! --- outputs
694  real (kind=kind_phys), dimension(:), intent(out) :: sfcemis
695 
696 ! --- locals:
697  integer :: i, i1, i2, j1, j2, idx
698 
699  real (kind=kind_phys) :: dltg, hdlt, tmp1, tmp2, &
700  & asnow, argh, hrgh, fsno, fsno0, fsno1
701 
702 ! --- reference emiss value for diff surface emiss index
703 ! 1-open water, 2-grass/shrub land, 3-bare soil, tundra,
704 ! 4-sandy desert, 5-rocky desert, 6-forest, 7-ice, 8-snow
705 
706  real (kind=kind_phys) :: emsref(8)
707  data emsref / 0.97, 0.95, 0.94, 0.90, 0.93, 0.96, 0.96, 0.99 /
708 
709 !
710 !===> ... begin here
711 !
713  if ( iemslw == 0 ) then ! sfc emiss default to 1.0
714 
715  sfcemis(:) = f_one
716  return
717 
718  else ! emiss set by sfc type and condition
719 
720  dltg = 360.0 / float(imxems)
721  hdlt = 0.5 * dltg
722 
723 ! --- ... mapping input data onto model grid
724 ! note: this is a simple mapping method, an upgrade is needed if
725 ! the model grid is much corcer than the 1-deg data resolution
726 
727  lab_do_imax : do i = 1, imax
728 
729  if ( nint(slmsk(i)) == 0 ) then ! sea point
730 
731  sfcemis(i) = emsref(1)
732 
733  else if ( nint(slmsk(i)) == 2 ) then ! sea-ice
734 
735  sfcemis(i) = emsref(7)
736 
737  else ! land
738 
739 ! --- map grid in longitude direction
740  i2 = 1
741  j2 = 1
742 
743  tmp1 = xlon(i) * rad2dg
744  if (tmp1 < f_zero) tmp1 = tmp1 + 360.0
745 
746  lab_do_imxems : do i1 = 1, imxems
747  tmp2 = dltg * (i1 - 1) + hdlt
748 
749  if (abs(tmp1-tmp2) <= hdlt) then
750  i2 = i1
751  exit lab_do_imxems
752  endif
753  enddo lab_do_imxems
754 
755 ! --- map grid in latitude direction
756  tmp1 = xlat(i) * rad2dg ! if xlat in pi/2 -> -pi/2 range
757 ! tmp1 = 90.0 - xlat(i)*rad2dg ! if xlat in 0 -> pi range
758 
759  lab_do_jmxems : do j1 = 1, jmxems
760  tmp2 = 90.0 - dltg * (j1 - 1)
761 
762  if (abs(tmp1-tmp2) <= hdlt) then
763  j2 = j1
764  exit lab_do_jmxems
765  endif
766  enddo lab_do_jmxems
767 
768 
769  idx = max( 2, idxems(i2,j2) )
770  if ( idx >= 7 ) idx = 2
771  sfcemis(i) = emsref(idx)
772 
773  endif ! end if_slmsk_block
774 
776 
777  if ( ialbflg==1 .and. nint(slmsk(i))==1 ) then ! input land area snow cover
778 
779  fsno0 = sncovr(i)
780  fsno1 = f_one - fsno0
781  sfcemis(i) = sfcemis(i)*fsno1 + emsref(8)*fsno0
782 
783  else ! compute snow cover from snow depth
784  if ( snowf(i) > f_zero ) then
785  asnow = 0.02*snowf(i)
786  argh = min(0.50, max(.025, 0.01*zorlf(i)))
787  hrgh = min(f_one, max(0.20, 1.0577-1.1538e-3*hprif(i) ) )
788  fsno0 = asnow / (argh + asnow) * hrgh
789  if (nint(slmsk(i)) == 0 .and. tsknf(i) > 271.2) &
790  & fsno0=f_zero
791  fsno1 = f_one - fsno0
792  sfcemis(i) = sfcemis(i)*fsno1 + emsref(8)*fsno0
793  endif
794 
795  endif ! end if_ialbflg
796 
797  enddo lab_do_imax
798 
799  endif ! end if_iemslw_block
800 
801 !chk print *,' In setemis, iemsflg, sfcemis =',iemsflg,sfcemis
802 
803 !
804  return
805 !...................................
806  end subroutine setemis
807 !-----------------------------------
808 
810 !
811 !.........................................!
812  end module module_radiation_surface !
813 !=========================================!
real(kind=kind_phys), parameter con_pi
pi
Definition: physcons.f:48
real(kind=kind_phys), parameter con_tice
temp freezing sea (K)
Definition: physcons.f:100
subroutine, public setemis(xlon, xlat, slmsk, snowf, sncovr, zorlf, tsknf, tairf, hprif, IMAX, sfcemis )
This subroutine computes surface emissivity for LW radiation.
integer, save ialbflg
surface albedo scheme control flag =0:vegetation type based climatological albedo scheme =1:seaso...
Definition: physparam.f:273
integer, dimension(:,:), allocatable idxems
global surface emissivity index array
subroutine, public sfc_init(me)
This subroutine is the initialization program for surface radiation related quantities (albedo...
character, save semis_file
external sfc emissivity data table: sfc_emissivity_idx.txt
Definition: physparam.f:281
integer, parameter, public nf_albd
num of sfc albedo components
real(kind=kind_phys), parameter con_t0c
temp at 0C (K)
Definition: physcons.f:96
integer, save iemsflg
surface emissivity scheme control flag =0:black-body surface emissivity(=1.0) =1:vegetation type ...
Definition: physparam.f:278
subroutine, public setalb(slmsk, snowf, sncovr, snoalb, zorlf, coszf, tsknf, tairf, hprif, alvsf, alnsf, alvwf, alnwf, facsf, facwf, fice, tisfc, IMAX, sfcalb )
This subroutine computes four components of surface albedos (i.e., vis-nir, direct-diffused) accordin...
integer iemslw
global surface emissivity contrl flag set up in 'sfc_init'
integer, parameter, public imxems
num of longitude points in global emis-type map
real(kind=kind_phys), parameter con_ttp
temp at H2O 3pt (K)
Definition: physcons.f:98
integer, parameter, public jmxems
num of latitude points in global emis-type map