CCPP SciDoc v7.0.x  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
radiation_surface.f
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,snodi,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! ( lsm,lsm_noahmp,lsm_ruc,frac_grid,cplice,use_lake_model, !
31! --- inputs:
32! lakefrac,xlon,xlat,slmsk,snodl,snodi,sncovr,sncovr_ice, !
33! zorlf,tsknf,tairf,hprif, !
34! semis_lnd,semis_ice,semis_wat,IMAX,fracl,fraco,fraci,icy, !
35!
36! --- outputs:
37! semisbase, sfcemis !
38!
39! !
40! external modules referenced: !
41! !
42! 'module machine' in 'machine.f' !
43! 'module physcons' in 'physcons.f' !
44! 'module module_iounitdef' in 'iounitdef.f' !
45! !
46! !
47! program history log: !
48! 1995 y.t. hou - created albaer.f (include albedo !
49! and aerosols calculations) !
50! nov 1997 y.t. hou - modified snow albedo !
51! jan 1998 y.t. hou - included grumbine's sea-ice scheme !
52! feb 1998 h.l. pan - seasonal interpolation in cycle !
53! mar 2000 y.t. hou - modified to use opac aerosol data !
54! apr 2003 y.t. hou - seperate albedo and aerosols into !
55! two subroutines, rewritten in f90 modulized form !
56! jan 2005 s. moorthi - xingren's sea-ice fraction added !
57! apr 2005 y.t. hou - revised module structure !
58! feb 2006 y.t. hou - add varying surface emissivity, !
59! modified sfc albedo structure for modis shceme !
60! Mar 2006 s. moorthi - added surface temp over ice fraction !
61! mar 2007 c. marshall & h. wei !
62! - added modis based sfc albedo scheme !
63! may 2007 y. hou & s. moorthi !
64! - fix bug in modis albedo over ocean !
65! aug 2007 h. wei & s. moorthi !
66! - fix bug in modis albedo over sea-ice !
67! aug 2007 y. hou - fix bug in emissivity over ocean in !
68! the modis scheme option !
69! dec 2008 f. yang - modified zenith angle dependence on !
70! surface albedo over land. (2008 jamc)!
71! aug 2012 y. hou - minor modification in initialization !
72! subr 'sfc_init'. !
73! nov 2012 y. hou - modified control parameters through !
74! module 'physparam'. !
75! jun 2018 h-m lin/y-t hou - correct error in clim-scheme of !
76! weak/strong factor and restore to the orig form !
77! !
78!!!!! ========================================================== !!!!!
79!!!!! end descriptions !!!!!
80!!!!! ========================================================== !!!!!
81
82
102
106!
107 use machine, only : kind_phys
108 use module_iounitdef, only : niradsf
109 use surface_perturbation, only : ppfbet
110!
111 implicit none
112!
113 private
114
115! --- version tag and last revision date
116 character(40), parameter :: &
117 & VTAGSFC='NCEP-Radiation_surface v5.1 Nov 2012 '
118! & VTAGSFC='NCEP-Radiation_surface v5.0 Aug 2012 '
119
120! --- constant parameters
121 integer, parameter, public :: imxems = 360 ! number of longtitude points in global emis-type map
122 integer, parameter, public :: jmxems = 180 ! number of latitude points in global emis-type map
123 real (kind=kind_phys), parameter :: f_zero = 0.0
124 real (kind=kind_phys), parameter :: f_one = 1.0
125 real (kind=kind_phys), parameter :: epsln = 1.0e-6
126 real (kind=kind_phys) :: rad2dg
127 integer, allocatable :: idxems(:,:) ! global surface emissivity index array
128 integer :: iemslw = 1 ! global surface emissivity control flag set up in 'sfc_init'
129!
130 public sfc_init, setalb, setemis
131 public f_zero, f_one, epsln
132
133! =================
134 contains
135! =================
136
140!-----------------------------------
141 subroutine sfc_init &
142 & ( me, ialbflg, iemsflg, semis_file, con_pi, errmsg, errflg )! --- inputs/outputs:
143!
144! =================================================================== !
145! !
146! this program is the initialization program for surface radiation !
147! related quantities (albedo, emissivity, etc.) !
148! !
149! usage: call sfc_init !
150! !
151! subprograms called: none !
152! !
153! ==================== defination of variables ==================== !
154! !
155! inputs: !
156! me - print control flag !
157! ialbflg - control flag for surface albedo schemes !
158! =1: use modis based surface albedo !
159! =2: use surface albedo from land model !
160! iemsflg - control flag for sfc emissivity schemes (ab:2-dig)!
161! a:=0 set sfc air/ground t same for lw radiation !
162! =1 set sfc air/ground t diff for lw radiation !
163! b:=1 use varying climtology sfc emiss (veg based) !
164! =2 use surface emissivity from land model !
165! !
166! outputs: (CCPP error handling) !
167! errmsg - CCPP error message !
168! errflg - CCPP error flag !
169! !
170! ==================== end of description ===================== !
171!
172 implicit none
173
174! --- inputs:
175 integer, intent(in) :: me, ialbflg, iemsflg
176 real(kind=kind_phys), intent(in) :: con_pi
177 character(len=26), intent(in) :: semis_file
178! --- outputs: ( none )
179 character(len=*), intent(out) :: errmsg
180 integer, intent(out) :: errflg
181
182! --- locals:
183 integer :: i, k
184! integer :: ia, ja
185 logical :: file_exist
186 character :: cline*80
187!
188!===> ... begin here
189!
190 errmsg = ''
191 errflg = 0
192!
193 ! Module
194 rad2dg = 180.0 / con_pi
195
196 if ( me == 0 ) print *, vtagsfc ! print out version tag
197
202
203 if ( ialbflg == 1 ) then
204
205 if ( me == 0 ) then
206 print *,' - Using MODIS based land surface albedo for sw'
207 endif
208
209 elseif ( ialbflg == 2 ) then ! use albedo from land model
210
211 if ( me == 0 ) then
212 print *,' - Using Albedo From Land Model'
213 endif
214
215 else
216
217 errmsg = 'module_radiation_surface: invalid ialbflg option'
218 errflg = 1
219 return
220
221 endif ! end if_ialbflg_block
222
227
228 iemslw = mod(iemsflg, 10) ! emissivity control
229
230 if ( iemslw == 1 ) then ! input sfc emiss type map
231
232! --- allocate data space
233 if ( .not. allocated(idxems) ) then
234 allocate ( idxems(imxems,jmxems) )
235 endif
236
237! --- check to see if requested emissivity data file existed
238
239 inquire (file=semis_file, exist=file_exist)
240
241 if ( .not. file_exist ) then
242 if ( me == 0 ) then
243 print *,' - Using Varying Surface Emissivity for lw'
244 print *,' Requested data file "',semis_file,'" not found!'
245 endif
246 errmsg = 'module_radiation_surface: surface emissivity
247 & file not provided'
248 errflg = 1
249 return
250
251 else
252 close(niradsf)
253 open (niradsf,file=semis_file,form='formatted',status='old')
254 rewind niradsf
255
256 read (niradsf,12) cline
257 12 format(a80)
258
259 read (niradsf,14) idxems
260 14 format(80i1)
261
262 if ( me == 0 ) then
263 print *,' - Using Varying Surface Emissivity for lw'
264 print *,' Opened data file: ',semis_file
265 print *, cline
266!check print *,' CHECK: Sample emissivity index data'
267! ia = IMXEMS / 5
268! ja = JMXEMS / 5
269! print *, idxems(1:IMXEMS:ia,1:JMXEMS:ja)
270 endif
271
272 close(niradsf)
273 endif ! end if_file_exist_block
274
275 elseif ( iemslw == 2 ) then ! use emiss from land model
276
277 if ( me == 0 ) then
278 print *,' - Using Surface Emissivity From Land Model'
279 endif
280
281 else
282
283 errmsg = 'module_radiation_surface: invalid iemslw option'
284 errflg = 1
285 return
286
287 endif ! end if_iemslw_block
288
289!
290 return
291!...................................
292 end subroutine sfc_init
293!-----------------------------------
294
343!-----------------------------------
344 subroutine setalb &
345 & ( slmsk,lsm,lsm_noahmp,lsm_ruc,use_cice_alb,snodi, & ! --- inputs:
346 & sncovr,sncovr_ice,snoalb,zorlf,coszf, &
347 & tsknf,tairf,hprif,frac_grid, lakefrac, &
348 & alvsf,alnsf,alvwf,alnwf,facsf,facwf,fice,tisfc, &
349 & lsmalbdvis, lsmalbdnir, lsmalbivis, lsmalbinir, &
350 & icealbdvis, icealbdnir, icealbivis, icealbinir, &
351 & imax, albppert, pertalb, fracl, fraco, fraci, icy, &
352 & ialbflg, con_ttp, &
353 & sfcalb & ! --- outputs:
354 & )
355
356! =================================================================== !
357! !
358! this program computes four components of surface albedos (i.e. !
359! vis-nir, direct-diffused) according to controflag ialbflg. !
360! 1) climatological surface albedo scheme (briegleb 1992) !
361! 2) modis retrieval based scheme from boston univ. !
362! !
363! !
364! usage: call setalb !
365! !
366! subprograms called: none !
367! !
368! ==================== defination of variables ==================== !
369! !
370! inputs: !
371! slmsk (IMAX) - sea(0),land(1),ice(2) mask on fcst model grid !
372! snodi (IMAX) - snow depth water equivalent in mm over ice !
373! sncovr(IMAX) - ialgflg=0: not used !
374! ialgflg=1: snow cover over land in fraction !
375! sncovr_ice(IMAX) - ialgflg=0: not used !
376! ialgflg=1: snow cover over ice in fraction !
377! snoalb(IMAX) - ialbflg=0: not used !
378! ialgflg=1: max snow albedo over land in fraction !
379! zorlf (IMAX) - surface roughness in cm !
380! coszf (IMAX) - cosin of solar zenith angle !
381! tsknf (IMAX) - ground surface temperature in k !
382! tairf (IMAX) - lowest model layer air temperature in k !
383! hprif (IMAX) - topographic sdv in m !
384! --- for ialbflg=0 climtological albedo scheme --- !
385! alvsf (IMAX) - 60 degree vis albedo with strong cosz dependency !
386! alnsf (IMAX) - 60 degree nir albedo with strong cosz dependency !
387! alvwf (IMAX) - 60 degree vis albedo with weak cosz dependency !
388! alnwf (IMAX) - 60 degree nir albedo with weak cosz dependency !
389! --- for ialbflg=1 modis based land albedo scheme --- !
390! alvsf (IMAX) - visible black sky albedo at zenith 60 degree !
391! alnsf (IMAX) - near-ir black sky albedo at zenith 60 degree !
392! alvwf (IMAX) - visible white sky albedo !
393! alnwf (IMAX) - near-ir white sky albedo !
394! !
395! facsf (IMAX) - fractional coverage with strong cosz dependency !
396! facwf (IMAX) - fractional coverage with weak cosz dependency !
397! fice (IMAX) - sea-ice fraction !
398! tisfc (IMAX) - sea-ice surface temperature !
399! IMAX - array horizontal dimension !
400! ialbflg - control flag for surface albedo schemes !
401! =1: use modis based surface albedo !
402! =2: use surface albedo from land model !
403! !
404! outputs: !
405! sfcalb(IMAX,NF_ALBD) !
406! ( :, 1) - near ir direct beam albedo !
407! ( :, 2) - near ir diffused albedo !
408! ( :, 3) - uv+vis direct beam albedo !
409! ( :, 4) - uv+vis diffused albedo !
410! !
411! ==================== end of description ===================== !
412!
413 implicit none
414
415! --- inputs
416 integer, intent(in) :: imax, ialbflg
417 integer, intent(in) :: lsm, lsm_noahmp, lsm_ruc
418 logical, intent(in) :: use_cice_alb, frac_grid
419
420 real (kind=kind_phys), dimension(:), intent(in) :: &
421 & lakefrac, &
422 & slmsk, snodi, zorlf, coszf, tsknf, tairf, hprif, &
423 & alvsf, alnsf, alvwf, alnwf, facsf, facwf, fice, tisfc, &
424 & sncovr, sncovr_ice, snoalb, albppert ! sfc-perts, mgehne
425 real (kind=kind_phys), dimension(:), intent(in), optional :: &
426 & icealbdvis, icealbdnir, icealbivis, icealbinir
427 real (kind=kind_phys), intent(in) :: pertalb, con_ttp! sfc-perts, mgehne
428 real (kind=kind_phys), dimension(:), intent(in) :: &
429 & fracl, fraco, fraci
430 real (kind=kind_phys), dimension(:),intent(inout) :: &
431 & lsmalbdvis, lsmalbdnir, lsmalbivis, lsmalbinir
432
433 logical, dimension(:), intent(in) :: &
434 & icy
435
436! --- outputs
437 real (kind=kind_phys), dimension(:,:), intent(out) :: sfcalb
438
439! --- locals:
440 real (kind=kind_phys) :: asnvb, asnnb, asnvd, asnnd, asevb &
441 &, asenb, asevd, asend, fsno, fsea, rfcs, rfcw, flnd &
442 &, asnow, argh, hrgh, fsno0, fsno1, flnd0, fsea0, csnow &
443 &, a1, a2, b1, b2, b3, ab1bm, ab2bm, m, s, alpha, beta, albtmp
444
445 real (kind=kind_phys) :: asevb_wat,asenb_wat,asevd_wat,asend_wat, &
446 & asevb_ice,asenb_ice,asevd_ice,asend_ice
447
448 real (kind=kind_phys) :: alndnb, alndnd, alndvb, alndvd
449
450 real (kind=kind_phys) ffw, dtgd, icealb
451 real (kind=kind_phys), parameter :: epsln=1.0e-8_kind_phys
452
453 integer :: i, k, kk, iflag
454
455!
456!===> ... begin here
457!
459 if ( ialbflg == 1 ) then
460
461 do i = 1, imax
462
463 !-- water albedo
464 asevd_wat = 0.06
465 asend_wat = 0.06
466 asevb_wat = asevd_wat
467 asenb_wat = asevd_wat
468
469 ! direct albedo CZA dependence over water
470 if (fraco(i) > f_zero .and. coszf(i) > 0.0001) then
471 asevb_wat = max(asevd_wat, 0.026/(coszf(i)**1.7 + 0.065) &
472 & + 0.15 * (coszf(i)-0.1) * (coszf(i)-0.5) &
473 & * (coszf(i)-f_one))
474 asenb_wat = asevb_wat
475 endif
476
477 if (icy(i)) then !-- Computation of ice albedo
478
479 if (use_cice_alb .and. lakefrac(i) < epsln) then
480 icealb = icealbivis(i)
481 else
482 icealb = f_zero
483 endif
484 if (icealb > epsln) then !-- use ice albedo from CICE for sea-ice
485 asevd_ice = icealbivis(i)
486 asend_ice = icealbinir(i)
487 asevb_ice = icealbdvis(i)
488 asenb_ice = icealbdnir(i)
489 else
490 asnow = 0.02*snodi(i)
491 argh = min(0.50, max(.025, 0.01*zorlf(i)))
492 hrgh = min(f_one,max(0.20,1.0577-1.1538e-3*hprif(i)))
493 fsno0 = asnow / (argh + asnow) * hrgh ! snow fraction on ice
494 ! diffused
495 if (tsknf(i) > 271.1 .and. tsknf(i) < 271.5) then
496 !tgs: looks like albedo reduction from puddles on ice
497 a1 = (tsknf(i) - 271.1)**2
498 asevd_ice = 0.7 - 4.0*a1
499 asend_ice = 0.65 - 3.6875*a1
500 else
501 asevd_ice = 0.70
502 asend_ice = 0.65
503 endif
504 ! direct
505 asevb_ice = asevd_ice
506 asenb_ice = asend_ice
507
508 if (fsno0 > f_zero) then ! Snow on ice
509 dtgd = max(f_zero, min(5.0, (con_ttp-tisfc(i)) ))
510 b1 = 0.03 * dtgd
511 asnvd = (asevd_ice + b1) ! diffused snow albedo
512 asnnd = (asend_ice + b1)
513 if (coszf(i) > 0.0001 .and. coszf(i) < 0.5) then ! direct snow albedo
514 csnow = 0.5 * (3.0 / (f_one+4.0*coszf(i)) - f_one)
515 asnvb = min( 0.98, asnvd+(f_one-asnvd)*csnow )
516 asnnb = min( 0.98, asnnd+(f_one-asnnd)*csnow )
517 else
518 asnvb = asnvd
519 asnnb = asnnd
520 endif
521
522 ! composite ice and snow albedos
523 asevd_ice = asevd_ice * (1. - fsno0) + asnvd * fsno0
524 asend_ice = asend_ice * (1. - fsno0) + asnnd * fsno0
525 asevb_ice = asevb_ice * (1. - fsno0) + asnvb * fsno0
526 asenb_ice = asenb_ice * (1. - fsno0) + asnnb * fsno0
527 endif ! snow
528 endif ! if (use_cice_alb .and. lakefrac < epsln)
529 else ! icy = false, fill in values
530 asevd_ice = 0.70
531 asend_ice = 0.65
532 asevb_ice = 0.70
533 asenb_ice = 0.65
534 endif ! end icy
535
536 if (fracl(i) > f_zero) then
539
540 fsno0 = sncovr(i) ! snow fraction on land
541
542 fsno1 = f_one - fsno0
543 flnd0 = min(f_one, facsf(i)+facwf(i))
544 flnd = flnd0 * fsno1 ! snow-free fraction
545 fsno = f_one - flnd ! snow-covered fraction
546
547 ! - use Fanglin's zenith angle treatment.
548 if (coszf(i) > 0.0001) then
549 rfcs = 1.775/(1.0+1.55*coszf(i))
550 else
551 !- no sun
552 rfcs = f_one
553 endif
554 !- zenith dependence is applied only to direct beam albedo
555 ab1bm = min(0.99, alnsf(i)*rfcs)
556 ab2bm = min(0.99, alvsf(i)*rfcs)
557
558 alndnb = ab1bm *flnd + snoalb(i) * fsno
559 alndnd = alnwf(i)*flnd + snoalb(i) * fsno
560 alndvb = ab2bm *flnd + snoalb(i) * fsno
561 alndvd = alvwf(i)*flnd + snoalb(i) * fsno
562 lsmalbdnir(i) = min(0.99,max(0.01,alndnb))
563 lsmalbinir(i) = min(0.99,max(0.01,alndnd))
564 lsmalbdvis(i) = min(0.99,max(0.01,alndvb))
565 lsmalbivis(i) = min(0.99,max(0.01,alndvd))
566 else
567 !-- fill in values for land albedo
568 alndnb = 0.
569 alndnd = 0.
570 alndvb = 0.
571 alndvd = 0.
572 endif ! end land
573
574 !-- Composite mean surface albedo from land, open water and
575 !-- ice fractions
576 sfcalb(i,1) = min(0.99,max(0.01,alndnb))*fracl(i) & ! direct beam NIR
577 & + asenb_wat*fraco(i) + asenb_ice*fraci(i)
578 sfcalb(i,2) = min(0.99,max(0.01,alndnd))*fracl(i) & ! diffuse NIR
579 & + asend_wat*fraco(i) + asend_ice*fraci(i)
580 sfcalb(i,3) = min(0.99,max(0.01,alndvb))*fracl(i) & ! direct beam visible
581 & + asevb_wat*fraco(i) + asevb_ice*fraci(i)
582 sfcalb(i,4) = min(0.99,max(0.01,alndvd))*fracl(i) & ! diffuse visible
583 & + asevd_wat*fraco(i) + asevd_ice*fraci(i)
584
585 enddo ! end_do_i_loop
586
588 elseif ( ialbflg == 2 ) then
589 do i = 1, imax
590
591 !-- water albedo
592 asevd_wat = 0.06
593 asend_wat = 0.06
594 asevb_wat = asevd_wat
595 asenb_wat = asevd_wat
596
597 ! direct albedo CZA dependence over water
598 if (fraco(i) > f_zero .and. coszf(i) > 0.0001) then
599 asevb_wat = max(asevd_wat, 0.026/(coszf(i)**1.7 + 0.065) &
600 & + 0.15 * (coszf(i)-0.1) * (coszf(i)-0.5) &
601 & * (coszf(i)-f_one))
602 asenb_wat = asevb_wat
603 endif
604
605 !-- ice albedo
606 !tgs: this part of the code needs the input from the ice
607 ! model. Otherwise it uses the backup albedo computation
608 ! from ialbflg = 1.
609
610 if (icy(i)) then !-- Computation of ice albedo
611
612 if (use_cice_alb .and. lakefrac(i) < epsln) then
613 icealb = icealbivis(i)
614 else
615 icealb = f_zero
616 endif
617
618 if (lsm == lsm_ruc .or. icealb > epsln) then !-- use ice albedo from the RUC ice model or
619 !-- use ice albedo from CICE for sea-ice
620 asevd_ice = icealbivis(i)
621 asend_ice = icealbinir(i)
622 asevb_ice = icealbdvis(i)
623 asenb_ice = icealbdnir(i)
624 else
625 !-- Computation of ice albedo
626 asnow = 0.02*snodi(i)
627 argh = min(0.50, max(.025, 0.01*zorlf(i)))
628 hrgh = min(f_one,max(0.20,1.0577-1.1538e-3*hprif(i)))
629 fsno0 = asnow / (argh + asnow) * hrgh
630 ! diffused
631 if (tsknf(i) > 271.1 .and. tsknf(i) < 271.5) then
632 !tgs: looks like albedo reduction from puddles on ice
633 a1 = (tsknf(i) - 271.1)**2
634 asevd_ice = 0.7 - 4.0*a1
635 asend_ice = 0.65 - 3.6875*a1
636 else
637 asevd_ice = 0.70
638 asend_ice = 0.65
639 endif
640 ! direct
641 asevb_ice = asevd_ice
642 asenb_ice = asend_ice
643
644 if (fsno0 > f_zero) then
645 ! Snow on ice
646 dtgd = max(f_zero, min(5.0, (con_ttp-tisfc(i)) ))
647 b1 = 0.03 * dtgd
648 asnvd = (asevd_ice + b1) ! diffused snow albedo
649 asnnd = (asend_ice + b1)
650
651 if (coszf(i) > 0.0001 .and. coszf(i) < 0.5) then ! direct snow albedo
652 csnow = 0.5 * (3.0 / (f_one+4.0*coszf(i)) - f_one)
653 asnvb = min( 0.98, asnvd+(f_one-asnvd)*csnow )
654 asnnb = min( 0.98, asnnd+(f_one-asnnd)*csnow )
655 else
656 asnvb = asnvd
657 asnnb = asnnd
658 endif
659
660 ! composite ice and snow albedos
661 asevd_ice = asevd_ice * (1. - fsno0) + asnvd * fsno0
662 asend_ice = asend_ice * (1. - fsno0) + asnnd * fsno0
663 asevb_ice = asevb_ice * (1. - fsno0) + asnvb * fsno0
664 asenb_ice = asenb_ice * (1. - fsno0) + asnnb * fsno0
665 endif ! snow
666 endif ! ice option from LSM or otherwise
667 else
668 ! icy = false, fill in values
669 asevd_ice = 0.70
670 asend_ice = 0.65
671 asevb_ice = 0.70
672 asenb_ice = 0.65
673 endif ! end icy
674
675 !-- Composite mean surface albedo from land, open water and
676 !-- ice fractions
677 sfcalb(i,1) = min(0.99,max(0.01,lsmalbdnir(i)))*fracl(i) & ! direct beam NIR
678 & + asenb_wat*fraco(i) + asenb_ice*fraci(i)
679 sfcalb(i,2) = min(0.99,max(0.01,lsmalbinir(i)))*fracl(i) & ! diffuse NIR
680 & + asend_wat*fraco(i) + asend_ice*fraci(i)
681 sfcalb(i,3) = min(0.99,max(0.01,lsmalbdvis(i)))*fracl(i) & ! direct beam visible
682 & + asevb_wat*fraco(i) + asevb_ice*fraci(i)
683 sfcalb(i,4) = min(0.99,max(0.01,lsmalbivis(i)))*fracl(i) & ! diffuse visible
684 & + asevd_wat*fraco(i) + asevd_ice*fraci(i)
685
686 enddo ! end_do_i_loop
687
688 endif ! end if_ialbflg
689!
690
691! sfc-perts, mgehne ***
693 if (pertalb>0.0) then
694 do i = 1, imax
695 do kk=1, 4
696 ! compute beta distribution parameters for all 4 albedos
697 m = sfcalb(i,kk)
698 s = pertalb*m*(1.-m)
699 alpha = m*m*(1.-m)/(s*s)-m
700 beta = alpha*(1.-m)/m
701 ! compute beta distribution value corresponding
702 ! to the given percentile albPpert to use as new albedo
703 call ppfbet(albppert(i),alpha,beta,iflag,albtmp)
704 sfcalb(i,kk) = albtmp
705 enddo
706 enddo ! end_do_i_loop
707 endif
708
709! *** sfc-perts, mgehne
710
711
712 return
713!...................................
714 end subroutine setalb
715!-----------------------------------
716
749!-----------------------------------
750 subroutine setemis &
751 & ( lsm,lsm_noahmp,lsm_ruc,frac_grid,cplice,use_lake_model, & ! --- inputs:
752 & lakefrac,xlon,xlat,slmsk,snodl,snodi,sncovr,sncovr_ice, &
753 & zorlf,tsknf,tairf,hprif, &
754 & semis_lnd,semis_ice,semis_wat,imax,fracl,fraco,fraci,icy, &
755 & semisbase, sfcemis & ! --- outputs:
756 & )
757
758! =================================================================== !
759! !
760! this program computes surface emissivity for lw radiation. !
761! !
762! usage: call setemis !
763! !
764! subprograms called: none !
765! !
766! ==================== defination of variables ==================== !
767! !
768! inputs: !
769! cplice - logical, ".true." when coupled to an ice model !
770! xlon (IMAX) - longitude in radiance, ok for both 0->2pi or !
771! -pi -> +pi ranges !
772! xlat (IMAX) - latitude in radiance, default to pi/2 -> -pi/2 !
773! range, otherwise see in-line comment !
774! slmsk (IMAX) - sea(0),land(1),ice(2) mask on fcst model grid !
775! snodl (IMAX) - snow depth water equivalent in mm over land !
776! snodi (IMAX) - snow depth water equivalent in mm over ice !
777! sncovr(IMAX) - ialbflg=1: snow cover over land in fraction !
778! sncovr_ice(IMAX) - snow cover over ice in fraction !
779! zorlf (IMAX) - surface roughness in cm !
780! tsknf (IMAX) - ground surface temperature in k !
781! tairf (IMAX) - lowest model layer air temperature in k !
782! hprif (IMAX) - topographic sdv in m !
783! IMAX - array horizontal dimension !
784! !
785! inputs/outputs: !
786! semis_lnd (IMAX) - land emissivity !
787! semis_ice (IMAX) - ice emissivity !
788! semis_wat (IMAX) - water emissivity !
789! !
790! outputs: !
791! sfcemis(IMAX) - surface emissivity !
792! !
793! ------------------------------------------------------------------- !
794! !
795! surface type definations: !
796! 1. open water 2. grass/wood/shrub land !
797! 3. tundra/bare soil 4. sandy desert !
798! 5. rocky desert 6. forest !
799! 7. ice 8. snow !
800! !
801! input index data lon from 0 towards east, lat from n to s !
802! !
803! ==================== end of description ===================== !
804!
805
806 implicit none
807
808! --- inputs
809 integer, intent(in) :: imax
810 integer, intent(in) :: lsm, lsm_noahmp, lsm_ruc
811 logical, intent(in) :: frac_grid, cplice
812 integer, dimension(:), intent(in) :: use_lake_model
813 real (kind=kind_phys), dimension(:), intent(in) :: lakefrac
814
815 real (kind=kind_phys), dimension(:), intent(in) :: &
816 & xlon,xlat, slmsk, snodl, snodi, sncovr, sncovr_ice, &
817 & zorlf, tsknf, tairf, hprif
818 real (kind=kind_phys), dimension(:), intent(in) :: &
819 & fracl, fraco, fraci
820 real (kind=kind_phys), dimension(:), intent(inout) :: &
821 & semis_lnd, semis_ice, semis_wat
822 logical, dimension(:), intent(in) :: &
823 & icy
824
825! --- outputs
826 real (kind=kind_phys), dimension(:), intent(out) :: semisbase
827 real (kind=kind_phys), dimension(:), intent(out) :: sfcemis
828
829! --- locals:
830 integer :: i, i1, i2, j1, j2, idx
831 integer :: ivgtyp
832
833 real (kind=kind_phys) :: dltg, hdlt, tmp1, tmp2, &
834 & asnow, argh, hrgh, fsno
835 real (kind=kind_phys) :: sfcemis_land, sfcemis_ice
836
837! --- reference emiss value for diff surface emiss index
838! 1-open water, 2-grass/shrub land, 3-bare soil, tundra,
839! 4-sandy desert, 5-rocky desert, 6-forest, 7-ice, 8-snow
840
841 real (kind=kind_phys) :: emsref(8)
842 data emsref / 0.97, 0.95, 0.94, 0.90, 0.93, 0.96, 0.96, 0.99 /
843
844!
845!===> ... begin here
846!
848
849 semis_wat = emsref(1)
850 if ( iemslw == 1 ) then
851
852 dltg = 360.0 / float(imxems)
853 hdlt = 0.5 * dltg
854
855! --- ... mapping input data onto model grid
856! note: this is a simple mapping method, an upgrade is needed if
857! the model grid is much coarser than the 1-deg data resolution
858
859 lab_do_imax : do i = 1, imax
860
861 if (.not. cplice .or. lakefrac(i) > f_zero) then
862 semis_ice(i) = emsref(7)
863 endif
864 if (fracl(i) < epsln) then ! no land
865 if ( abs(fraco(i)-f_one) < epsln ) then ! open water point
866 sfcemis(i) = emsref(1)
867 elseif ( abs(fraci(i)-f_one) < epsln ) then ! complete sea/lake ice
868 sfcemis(i) = semis_ice(i)
869 else
870 !-- fractional sea ice
871 sfcemis(i) = fraco(i)*emsref(1) + fraci(i)*semis_ice(i)
872 endif
873
874 else ! land or fractional grid
875
876! --- map grid in longitude direction
877 i2 = 1
878 j2 = 1
879
880 tmp1 = xlon(i) * rad2dg
881 if (tmp1 < f_zero) tmp1 = tmp1 + 360.0
882
883 lab_do_imxems : do i1 = 1, imxems
884 tmp2 = dltg * (i1 - 1) + hdlt
885
886 if (abs(tmp1-tmp2) <= hdlt) then
887 i2 = i1
888 exit lab_do_imxems
889 endif
890 enddo lab_do_imxems
891
892! --- map grid in latitude direction
893 tmp1 = xlat(i) * rad2dg ! if xlat in pi/2 -> -pi/2 range
894! tmp1 = 90.0 - xlat(i)*rad2dg ! if xlat in 0 -> pi range
895
896 lab_do_jmxems : do j1 = 1, jmxems
897 tmp2 = 90.0 - dltg * (j1 - 1)
898
899 if (abs(tmp1-tmp2) <= hdlt) then
900 j2 = j1
901 exit lab_do_jmxems
902 endif
903 enddo lab_do_jmxems
904
905 idx = max( 2, idxems(i2,j2) )
906 if ( idx >= 7 ) idx = 2
907 if (abs(fracl(i)-f_one) < epsln) then
908 sfcemis(i) = emsref(idx)
909 else
910 sfcemis(i) = fracl(i)*emsref(idx) + fraco(i)*emsref(1) &
911 & + fraci(i)*emsref(7)
912 endif
913 semisbase(i) = sfcemis(i)
914 semis_lnd(i) = emsref(idx)
915
916 endif
917
921
922 if (fracl(i) > epsln) then
923 if (sncovr(i) > f_zero) then
924 semis_lnd(i) = semis_lnd(i) * (f_one - sncovr(i)) &
925 & + emsref(8) * sncovr(i)
926 elseif (snodl(i) > f_zero) then
927 asnow = 0.02*snodl(i)
928 argh = min(0.50, max(.025, 0.01*zorlf(i)))
929 hrgh = min(f_one, max(0.20, 1.0577-1.1538e-3*hprif(i) ) )
930 fsno = min(f_one, max(f_zero, asnow/(argh+asnow) * hrgh))
931 semis_lnd(i) = semis_lnd(i)*(f_one-fsno) + emsref(8)*fsno
932 endif
933 endif
934 if (fraci(i) > epsln .and. &
935 & (lakefrac(i) > f_zero .or. .not. cplice)) then
936 if (sncovr_ice(i) > f_zero) then
937 semis_ice(i) = semis_ice(i) * (f_one - sncovr_ice(i)) &
938 & + emsref(8) * sncovr_ice(i)
939 elseif (snodi(i) > f_zero) then
940 asnow = 0.02*snodi(i)
941 argh = min(0.50, max(.025, 0.01*zorlf(i)))
942 hrgh = min(f_one, max(0.20, 1.0577-1.1538e-3*hprif(i) ) )
943 fsno = min(f_one, max(f_zero, asnow/(argh+asnow) * hrgh))
944 semis_ice(i) = semis_ice(i)*(f_one-fsno) + emsref(8)*fsno
945 endif
946 endif
947 sfcemis(i) = fracl(i)*semis_lnd(i) + fraco(i)*emsref(1) &
948 & + fraci(i)*semis_ice(i)
949
950 enddo lab_do_imax
951
952 elseif ( iemslw == 2 ) then ! sfc emiss updated in land model: Noah MP or RUC
953
954 do i = 1, imax
955
956 sfcemis_ice = emsref(7)
957 if ( icy(i) ) then !-- ice emissivity
958
959 !-- complete or fractional ice
960 if (lsm == lsm_noahmp) then
961 if (.not. cplice .or. lakefrac(i) > f_zero) then
962 if (sncovr_ice(i) > f_zero) then
963 sfcemis_ice = emsref(7) * (f_one-sncovr_ice(i)) &
964 & + emsref(8) * sncovr_ice(i)
965 elseif (snodi(i) > f_zero) then
966 asnow = 0.02*snodi(i)
967 argh = min(0.50, max(.025,0.01*zorlf(i)))
968 hrgh = min(f_one,max(0.20,1.0577-1.1538e-3*hprif(i)))
969 fsno = asnow / (argh + asnow) * hrgh
970 sfcemis_ice = emsref(7)*(f_one-fsno) + emsref(8)*fsno
971 endif
972 semis_ice(i) = sfcemis_ice
973 else
974 sfcemis_ice = semis_ice(i) ! output from CICE
975 endif
976 elseif (lsm == lsm_ruc) then
977 if (use_lake_model(i)>0) then
978 if (sncovr_ice(i) > f_zero) then
979 sfcemis_ice = emsref(7) * (f_one-sncovr_ice(i)) &
980 & + emsref(8) * sncovr_ice(i)
981 elseif (snodi(i) > f_zero) then
982 asnow = 0.02*snodi(i)
983 argh = min(0.50, max(.025,0.01*zorlf(i)))
984 hrgh = min(f_one,max(0.20,1.0577-1.1538e-3*hprif(i)))
985 fsno = asnow / (argh + asnow) * hrgh
986 sfcemis_ice = emsref(7)*(f_one-fsno) + emsref(8)*fsno
987 endif
988 semis_ice(i) = sfcemis_ice
989 else
990 sfcemis_ice = semis_ice(i) ! output from CICE or from RUC lsm (with snow effect)
991 endif
992 endif ! lsm check
993 endif ! icy
994
995 !-- land emissivity
996 !-- from Noah MP or RUC lsms
997 sfcemis_land = semis_lnd(i) ! albedo with snow effect from LSM
998
999 !-- Composite emissivity from land, water and ice fractions.
1000 sfcemis(i) = fracl(i)*sfcemis_land + fraco(i)*emsref(1) &
1001 & + fraci(i)*sfcemis_ice
1002
1003 enddo ! i
1004
1005 endif ! end if_iemslw_block
1006
1007!chk print *,' In setemis, iemsflg, sfcemis =',iemsflg,sfcemis
1008
1009!
1010 return
1011!...................................
1012 end subroutine setemis
1013!-----------------------------------
1014
1015!.........................................!
1016 end module module_radiation_surface !
1018!=========================================!
subroutine csnow
This subroutine calculates snow termal conductivity.
Definition sflx.f:1229
subroutine, public ppfbet(pr, p, q, iflag, x)
This subroutine computes the beta distribution value that matches the percentile from the random patt...
real(kind=kind_phys), parameter, public f_zero
real(kind=kind_phys), parameter, public f_one
integer, parameter, public jmxems
subroutine, public sfc_init(me, ialbflg, iemsflg, semis_file, con_pi, errmsg, errflg)
This subroutine is the initialization program for surface radiation related quantities (albedo,...
subroutine, public setemis(lsm, lsm_noahmp, lsm_ruc, frac_grid, cplice, use_lake_model, lakefrac, xlon, xlat, slmsk, snodl, snodi, sncovr, sncovr_ice, zorlf, tsknf, tairf, hprif, semis_lnd, semis_ice, semis_wat, imax, fracl, fraco, fraci, icy, semisbase, sfcemis)
This subroutine computes surface emissivity for LW radiation.
subroutine, public setalb(slmsk, lsm, lsm_noahmp, lsm_ruc, use_cice_alb, snodi, sncovr, sncovr_ice, snoalb, zorlf, coszf, tsknf, tairf, hprif, frac_grid, lakefrac, alvsf, alnsf, alvwf, alnwf, facsf, facwf, fice, tisfc, lsmalbdvis, lsmalbdnir, lsmalbivis, lsmalbinir, icealbdvis, icealbdnir, icealbivis, icealbinir, imax, albppert, pertalb, fracl, fraco, fraci, icy, ialbflg, con_ttp, sfcalb)
This subroutine computes four components of surface albedos (i.e., vis-nir, direct-diffused) accordin...
integer, parameter, public imxems
integer, dimension(:,:), allocatable idxems
real(kind=kind_phys) rad2dg
real(kind=kind_phys), parameter, public epsln
This module contain the RUC land surface model driver.
Definition lsm_ruc.F90:5
this module defines fortran unit numbers for input/output data files for the ncep gfs model.
Definition iounitdef.f:53
This module sets up surface albedo for SW radiation and surface emissivity for LW radiation.
This module contains routines used in the percentile matching algorithm for the albedo and vegetation...