CCPP SciDoc v7.0.0  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, nf_albd, 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, nf_albd, 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(IMAX,NF_ALBD), intent(out) :: &
438 & sfcalb
439
440! --- locals:
441 real (kind=kind_phys) :: asnvb, asnnb, asnvd, asnnd, asevb &
442 &, asenb, asevd, asend, fsno, fsea, rfcs, rfcw, flnd &
443 &, asnow, argh, hrgh, fsno0, fsno1, flnd0, fsea0, csnow &
444 &, a1, a2, b1, b2, b3, ab1bm, ab2bm, m, s, alpha, beta, albtmp
445
446 real (kind=kind_phys) :: asevb_wat,asenb_wat,asevd_wat,asend_wat, &
447 & asevb_ice,asenb_ice,asevd_ice,asend_ice
448
449 real (kind=kind_phys) :: alndnb, alndnd, alndvb, alndvd
450
451 real (kind=kind_phys) ffw, dtgd, icealb
452 real (kind=kind_phys), parameter :: epsln=1.0e-8_kind_phys
453
454 integer :: i, k, kk, iflag
455
456!
457!===> ... begin here
458!
460 if ( ialbflg == 1 ) then
461
462 do i = 1, imax
463
464 !-- water albedo
465 asevd_wat = 0.06
466 asend_wat = 0.06
467 asevb_wat = asevd_wat
468 asenb_wat = asevd_wat
469
470 ! direct albedo CZA dependence over water
471 if (fraco(i) > f_zero .and. coszf(i) > 0.0001) then
472 asevb_wat = max(asevd_wat, 0.026/(coszf(i)**1.7 + 0.065) &
473 & + 0.15 * (coszf(i)-0.1) * (coszf(i)-0.5) &
474 & * (coszf(i)-f_one))
475 asenb_wat = asevb_wat
476 endif
477
478 if (icy(i)) then !-- Computation of ice albedo
479
480 if (use_cice_alb .and. lakefrac(i) < epsln) then
481 icealb = icealbivis(i)
482 else
483 icealb = f_zero
484 endif
485 if (icealb > epsln) then !-- use ice albedo from CICE for sea-ice
486 asevd_ice = icealbivis(i)
487 asend_ice = icealbinir(i)
488 asevb_ice = icealbdvis(i)
489 asenb_ice = icealbdnir(i)
490 else
491 asnow = 0.02*snodi(i)
492 argh = min(0.50, max(.025, 0.01*zorlf(i)))
493 hrgh = min(f_one,max(0.20,1.0577-1.1538e-3*hprif(i)))
494 fsno0 = asnow / (argh + asnow) * hrgh ! snow fraction on ice
495 ! diffused
496 if (tsknf(i) > 271.1 .and. tsknf(i) < 271.5) then
497 !tgs: looks like albedo reduction from puddles on ice
498 a1 = (tsknf(i) - 271.1)**2
499 asevd_ice = 0.7 - 4.0*a1
500 asend_ice = 0.65 - 3.6875*a1
501 else
502 asevd_ice = 0.70
503 asend_ice = 0.65
504 endif
505 ! direct
506 asevb_ice = asevd_ice
507 asenb_ice = asend_ice
508
509 if (fsno0 > f_zero) then ! Snow on ice
510 dtgd = max(f_zero, min(5.0, (con_ttp-tisfc(i)) ))
511 b1 = 0.03 * dtgd
512 asnvd = (asevd_ice + b1) ! diffused snow albedo
513 asnnd = (asend_ice + b1)
514 if (coszf(i) > 0.0001 .and. coszf(i) < 0.5) then ! direct snow albedo
515 csnow = 0.5 * (3.0 / (f_one+4.0*coszf(i)) - f_one)
516 asnvb = min( 0.98, asnvd+(f_one-asnvd)*csnow )
517 asnnb = min( 0.98, asnnd+(f_one-asnnd)*csnow )
518 else
519 asnvb = asnvd
520 asnnb = asnnd
521 endif
522
523 ! composite ice and snow albedos
524 asevd_ice = asevd_ice * (1. - fsno0) + asnvd * fsno0
525 asend_ice = asend_ice * (1. - fsno0) + asnnd * fsno0
526 asevb_ice = asevb_ice * (1. - fsno0) + asnvb * fsno0
527 asenb_ice = asenb_ice * (1. - fsno0) + asnnb * fsno0
528 endif ! snow
529 endif ! if (use_cice_alb .and. lakefrac < epsln)
530 else ! icy = false, fill in values
531 asevd_ice = 0.70
532 asend_ice = 0.65
533 asevb_ice = 0.70
534 asenb_ice = 0.65
535 endif ! end icy
536
537 if (fracl(i) > f_zero) then
540
541 fsno0 = sncovr(i) ! snow fraction on land
542
543 fsno1 = f_one - fsno0
544 flnd0 = min(f_one, facsf(i)+facwf(i))
545 flnd = flnd0 * fsno1 ! snow-free fraction
546 fsno = f_one - flnd ! snow-covered fraction
547
548 ! - use Fanglin's zenith angle treatment.
549 if (coszf(i) > 0.0001) then
550 rfcs = 1.775/(1.0+1.55*coszf(i))
551 else
552 !- no sun
553 rfcs = f_one
554 endif
555 !- zenith dependence is applied only to direct beam albedo
556 ab1bm = min(0.99, alnsf(i)*rfcs)
557 ab2bm = min(0.99, alvsf(i)*rfcs)
558
559 alndnb = ab1bm *flnd + snoalb(i) * fsno
560 alndnd = alnwf(i)*flnd + snoalb(i) * fsno
561 alndvb = ab2bm *flnd + snoalb(i) * fsno
562 alndvd = alvwf(i)*flnd + snoalb(i) * fsno
563 lsmalbdnir(i) = min(0.99,max(0.01,alndnb))
564 lsmalbinir(i) = min(0.99,max(0.01,alndnd))
565 lsmalbdvis(i) = min(0.99,max(0.01,alndvb))
566 lsmalbivis(i) = min(0.99,max(0.01,alndvd))
567 else
568 !-- fill in values for land albedo
569 alndnb = 0.
570 alndnd = 0.
571 alndvb = 0.
572 alndvd = 0.
573 endif ! end land
574
575 !-- Composite mean surface albedo from land, open water and
576 !-- ice fractions
577 sfcalb(i,1) = min(0.99,max(0.01,alndnb))*fracl(i) & ! direct beam NIR
578 & + asenb_wat*fraco(i) + asenb_ice*fraci(i)
579 sfcalb(i,2) = min(0.99,max(0.01,alndnd))*fracl(i) & ! diffuse NIR
580 & + asend_wat*fraco(i) + asend_ice*fraci(i)
581 sfcalb(i,3) = min(0.99,max(0.01,alndvb))*fracl(i) & ! direct beam visible
582 & + asevb_wat*fraco(i) + asevb_ice*fraci(i)
583 sfcalb(i,4) = min(0.99,max(0.01,alndvd))*fracl(i) & ! diffuse visible
584 & + asevd_wat*fraco(i) + asevd_ice*fraci(i)
585
586 enddo ! end_do_i_loop
587
589 elseif ( ialbflg == 2 ) then
590 do i = 1, imax
591
592 !-- water albedo
593 asevd_wat = 0.06
594 asend_wat = 0.06
595 asevb_wat = asevd_wat
596 asenb_wat = asevd_wat
597
598 ! direct albedo CZA dependence over water
599 if (fraco(i) > f_zero .and. coszf(i) > 0.0001) then
600 asevb_wat = max(asevd_wat, 0.026/(coszf(i)**1.7 + 0.065) &
601 & + 0.15 * (coszf(i)-0.1) * (coszf(i)-0.5) &
602 & * (coszf(i)-f_one))
603 asenb_wat = asevb_wat
604 endif
605
606 !-- ice albedo
607 !tgs: this part of the code needs the input from the ice
608 ! model. Otherwise it uses the backup albedo computation
609 ! from ialbflg = 1.
610
611 if (icy(i)) then !-- Computation of ice albedo
612
613 if (use_cice_alb .and. lakefrac(i) < epsln) then
614 icealb = icealbivis(i)
615 else
616 icealb = f_zero
617 endif
618
619 if (lsm == lsm_ruc .or. icealb > epsln) then !-- use ice albedo from the RUC ice model or
620 !-- use ice albedo from CICE for sea-ice
621 asevd_ice = icealbivis(i)
622 asend_ice = icealbinir(i)
623 asevb_ice = icealbdvis(i)
624 asenb_ice = icealbdnir(i)
625 else
626 !-- Computation of ice albedo
627 asnow = 0.02*snodi(i)
628 argh = min(0.50, max(.025, 0.01*zorlf(i)))
629 hrgh = min(f_one,max(0.20,1.0577-1.1538e-3*hprif(i)))
630 fsno0 = asnow / (argh + asnow) * hrgh
631 ! diffused
632 if (tsknf(i) > 271.1 .and. tsknf(i) < 271.5) then
633 !tgs: looks like albedo reduction from puddles on ice
634 a1 = (tsknf(i) - 271.1)**2
635 asevd_ice = 0.7 - 4.0*a1
636 asend_ice = 0.65 - 3.6875*a1
637 else
638 asevd_ice = 0.70
639 asend_ice = 0.65
640 endif
641 ! direct
642 asevb_ice = asevd_ice
643 asenb_ice = asend_ice
644
645 if (fsno0 > f_zero) then
646 ! Snow on ice
647 dtgd = max(f_zero, min(5.0, (con_ttp-tisfc(i)) ))
648 b1 = 0.03 * dtgd
649 asnvd = (asevd_ice + b1) ! diffused snow albedo
650 asnnd = (asend_ice + b1)
651
652 if (coszf(i) > 0.0001 .and. coszf(i) < 0.5) then ! direct snow albedo
653 csnow = 0.5 * (3.0 / (f_one+4.0*coszf(i)) - f_one)
654 asnvb = min( 0.98, asnvd+(f_one-asnvd)*csnow )
655 asnnb = min( 0.98, asnnd+(f_one-asnnd)*csnow )
656 else
657 asnvb = asnvd
658 asnnb = asnnd
659 endif
660
661 ! composite ice and snow albedos
662 asevd_ice = asevd_ice * (1. - fsno0) + asnvd * fsno0
663 asend_ice = asend_ice * (1. - fsno0) + asnnd * fsno0
664 asevb_ice = asevb_ice * (1. - fsno0) + asnvb * fsno0
665 asenb_ice = asenb_ice * (1. - fsno0) + asnnb * fsno0
666 endif ! snow
667 endif ! ice option from LSM or otherwise
668 else
669 ! icy = false, fill in values
670 asevd_ice = 0.70
671 asend_ice = 0.65
672 asevb_ice = 0.70
673 asenb_ice = 0.65
674 endif ! end icy
675
676 !-- Composite mean surface albedo from land, open water and
677 !-- ice fractions
678 sfcalb(i,1) = min(0.99,max(0.01,lsmalbdnir(i)))*fracl(i) & ! direct beam NIR
679 & + asenb_wat*fraco(i) + asenb_ice*fraci(i)
680 sfcalb(i,2) = min(0.99,max(0.01,lsmalbinir(i)))*fracl(i) & ! diffuse NIR
681 & + asend_wat*fraco(i) + asend_ice*fraci(i)
682 sfcalb(i,3) = min(0.99,max(0.01,lsmalbdvis(i)))*fracl(i) & ! direct beam visible
683 & + asevb_wat*fraco(i) + asevb_ice*fraci(i)
684 sfcalb(i,4) = min(0.99,max(0.01,lsmalbivis(i)))*fracl(i) & ! diffuse visible
685 & + asevd_wat*fraco(i) + asevd_ice*fraci(i)
686
687 enddo ! end_do_i_loop
688
689 endif ! end if_ialbflg
690!
691
692! sfc-perts, mgehne ***
694 if (pertalb>0.0) then
695 do i = 1, imax
696 do kk=1, 4
697 ! compute beta distribution parameters for all 4 albedos
698 m = sfcalb(i,kk)
699 s = pertalb*m*(1.-m)
700 alpha = m*m*(1.-m)/(s*s)-m
701 beta = alpha*(1.-m)/m
702 ! compute beta distribution value corresponding
703 ! to the given percentile albPpert to use as new albedo
704 call ppfbet(albppert(i),alpha,beta,iflag,albtmp)
705 sfcalb(i,kk) = albtmp
706 enddo
707 enddo ! end_do_i_loop
708 endif
709
710! *** sfc-perts, mgehne
711
712
713 return
714!...................................
715 end subroutine setalb
716!-----------------------------------
717
750!-----------------------------------
751 subroutine setemis &
752 & ( lsm,lsm_noahmp,lsm_ruc,frac_grid,cplice,use_lake_model, & ! --- inputs:
753 & lakefrac,xlon,xlat,slmsk,snodl,snodi,sncovr,sncovr_ice, &
754 & zorlf,tsknf,tairf,hprif, &
755 & semis_lnd,semis_ice,semis_wat,imax,fracl,fraco,fraci,icy, &
756 & semisbase, sfcemis & ! --- outputs:
757 & )
758
759! =================================================================== !
760! !
761! this program computes surface emissivity for lw radiation. !
762! !
763! usage: call setemis !
764! !
765! subprograms called: none !
766! !
767! ==================== defination of variables ==================== !
768! !
769! inputs: !
770! cplice - logical, ".true." when coupled to an ice model !
771! xlon (IMAX) - longitude in radiance, ok for both 0->2pi or !
772! -pi -> +pi ranges !
773! xlat (IMAX) - latitude in radiance, default to pi/2 -> -pi/2 !
774! range, otherwise see in-line comment !
775! slmsk (IMAX) - sea(0),land(1),ice(2) mask on fcst model grid !
776! snodl (IMAX) - snow depth water equivalent in mm over land !
777! snodi (IMAX) - snow depth water equivalent in mm over ice !
778! sncovr(IMAX) - ialbflg=1: snow cover over land in fraction !
779! sncovr_ice(IMAX) - snow cover over ice in fraction !
780! zorlf (IMAX) - surface roughness in cm !
781! tsknf (IMAX) - ground surface temperature in k !
782! tairf (IMAX) - lowest model layer air temperature in k !
783! hprif (IMAX) - topographic sdv in m !
784! IMAX - array horizontal dimension !
785! !
786! inputs/outputs: !
787! semis_lnd (IMAX) - land emissivity !
788! semis_ice (IMAX) - ice emissivity !
789! semis_wat (IMAX) - water emissivity !
790! !
791! outputs: !
792! sfcemis(IMAX) - surface emissivity !
793! !
794! ------------------------------------------------------------------- !
795! !
796! surface type definations: !
797! 1. open water 2. grass/wood/shrub land !
798! 3. tundra/bare soil 4. sandy desert !
799! 5. rocky desert 6. forest !
800! 7. ice 8. snow !
801! !
802! input index data lon from 0 towards east, lat from n to s !
803! !
804! ==================== end of description ===================== !
805!
808
809 implicit none
810
811! --- inputs
812 integer, intent(in) :: imax
813 integer, intent(in) :: lsm, lsm_noahmp, lsm_ruc
814 logical, intent(in) :: frac_grid, cplice
815 integer, dimension(:), intent(in) :: use_lake_model
816 real (kind=kind_phys), dimension(:), intent(in) :: lakefrac
817
818 real (kind=kind_phys), dimension(:), intent(in) :: &
819 & xlon,xlat, slmsk, snodl, snodi, sncovr, sncovr_ice, &
820 & zorlf, tsknf, tairf, hprif
821 real (kind=kind_phys), dimension(:), intent(in) :: &
822 & fracl, fraco, fraci
823 real (kind=kind_phys), dimension(:), intent(inout) :: &
824 & semis_lnd, semis_ice, semis_wat
825 logical, dimension(:), intent(in) :: &
826 & icy
827
828! --- outputs
829 real (kind=kind_phys), dimension(:), intent(out) :: semisbase
830 real (kind=kind_phys), dimension(:), intent(out) :: sfcemis
831
832! --- locals:
833 integer :: i, i1, i2, j1, j2, idx
834 integer :: ivgtyp
835
836 real (kind=kind_phys) :: dltg, hdlt, tmp1, tmp2, &
837 & asnow, argh, hrgh, fsno
838 real (kind=kind_phys) :: sfcemis_land, sfcemis_ice
839
840! --- reference emiss value for diff surface emiss index
841! 1-open water, 2-grass/shrub land, 3-bare soil, tundra,
842! 4-sandy desert, 5-rocky desert, 6-forest, 7-ice, 8-snow
843
844 real (kind=kind_phys) :: emsref(8)
845 data emsref / 0.97, 0.95, 0.94, 0.90, 0.93, 0.96, 0.96, 0.99 /
846
847!
848!===> ... begin here
849!
851
852 semis_wat = emsref(1)
853 if ( iemslw == 1 ) then
854
855 dltg = 360.0 / float(imxems)
856 hdlt = 0.5 * dltg
857
858! --- ... mapping input data onto model grid
859! note: this is a simple mapping method, an upgrade is needed if
860! the model grid is much coarser than the 1-deg data resolution
861
862 lab_do_imax : do i = 1, imax
863
864 if (.not. cplice .or. lakefrac(i) > f_zero) then
865 semis_ice(i) = emsref(7)
866 endif
867 if (fracl(i) < epsln) then ! no land
868 if ( abs(fraco(i)-f_one) < epsln ) then ! open water point
869 sfcemis(i) = emsref(1)
870 elseif ( abs(fraci(i)-f_one) < epsln ) then ! complete sea/lake ice
871 sfcemis(i) = semis_ice(i)
872 else
873 !-- fractional sea ice
874 sfcemis(i) = fraco(i)*emsref(1) + fraci(i)*semis_ice(i)
875 endif
876
877 else ! land or fractional grid
878
879! --- map grid in longitude direction
880 i2 = 1
881 j2 = 1
882
883 tmp1 = xlon(i) * rad2dg
884 if (tmp1 < f_zero) tmp1 = tmp1 + 360.0
885
886 lab_do_imxems : do i1 = 1, imxems
887 tmp2 = dltg * (i1 - 1) + hdlt
888
889 if (abs(tmp1-tmp2) <= hdlt) then
890 i2 = i1
891 exit lab_do_imxems
892 endif
893 enddo lab_do_imxems
894
895! --- map grid in latitude direction
896 tmp1 = xlat(i) * rad2dg ! if xlat in pi/2 -> -pi/2 range
897! tmp1 = 90.0 - xlat(i)*rad2dg ! if xlat in 0 -> pi range
898
899 lab_do_jmxems : do j1 = 1, jmxems
900 tmp2 = 90.0 - dltg * (j1 - 1)
901
902 if (abs(tmp1-tmp2) <= hdlt) then
903 j2 = j1
904 exit lab_do_jmxems
905 endif
906 enddo lab_do_jmxems
907
908 idx = max( 2, idxems(i2,j2) )
909 if ( idx >= 7 ) idx = 2
910 if (abs(fracl(i)-f_one) < epsln) then
911 sfcemis(i) = emsref(idx)
912 else
913 sfcemis(i) = fracl(i)*emsref(idx) + fraco(i)*emsref(1) &
914 & + fraci(i)*emsref(7)
915 endif
916 semisbase(i) = sfcemis(i)
917 semis_lnd(i) = emsref(idx)
918
919 endif
920
924
925 if (fracl(i) > epsln) then
926 if (sncovr(i) > f_zero) then
927 semis_lnd(i) = semis_lnd(i) * (f_one - sncovr(i)) &
928 & + emsref(8) * sncovr(i)
929 elseif (snodl(i) > f_zero) then
930 asnow = 0.02*snodl(i)
931 argh = min(0.50, max(.025, 0.01*zorlf(i)))
932 hrgh = min(f_one, max(0.20, 1.0577-1.1538e-3*hprif(i) ) )
933 fsno = min(f_one, max(f_zero, asnow/(argh+asnow) * hrgh))
934 semis_lnd(i) = semis_lnd(i)*(f_one-fsno) + emsref(8)*fsno
935 endif
936 endif
937 if (fraci(i) > epsln .and. &
938 & (lakefrac(i) > f_zero .or. .not. cplice)) then
939 if (sncovr_ice(i) > f_zero) then
940 semis_ice(i) = semis_ice(i) * (f_one - sncovr_ice(i)) &
941 & + emsref(8) * sncovr_ice(i)
942 elseif (snodi(i) > f_zero) then
943 asnow = 0.02*snodi(i)
944 argh = min(0.50, max(.025, 0.01*zorlf(i)))
945 hrgh = min(f_one, max(0.20, 1.0577-1.1538e-3*hprif(i) ) )
946 fsno = min(f_one, max(f_zero, asnow/(argh+asnow) * hrgh))
947 semis_ice(i) = semis_ice(i)*(f_one-fsno) + emsref(8)*fsno
948 endif
949 endif
950 sfcemis(i) = fracl(i)*semis_lnd(i) + fraco(i)*emsref(1) &
951 & + fraci(i)*semis_ice(i)
952
953 enddo lab_do_imax
954
955 elseif ( iemslw == 2 ) then ! sfc emiss updated in land model: Noah MP or RUC
956
957 do i = 1, imax
958
959 sfcemis_ice = emsref(7)
960 if ( icy(i) ) then !-- ice emissivity
961
962 !-- complete or fractional ice
963 if (lsm == lsm_noahmp) then
964 if (.not. cplice .or. lakefrac(i) > f_zero) then
965 if (sncovr_ice(i) > f_zero) then
966 sfcemis_ice = emsref(7) * (f_one-sncovr_ice(i)) &
967 & + emsref(8) * sncovr_ice(i)
968 elseif (snodi(i) > f_zero) then
969 asnow = 0.02*snodi(i)
970 argh = min(0.50, max(.025,0.01*zorlf(i)))
971 hrgh = min(f_one,max(0.20,1.0577-1.1538e-3*hprif(i)))
972 fsno = asnow / (argh + asnow) * hrgh
973 sfcemis_ice = emsref(7)*(f_one-fsno) + emsref(8)*fsno
974 endif
975 semis_ice(i) = sfcemis_ice
976 else
977 sfcemis_ice = semis_ice(i) ! output from CICE
978 endif
979 elseif (lsm == lsm_ruc) then
980 if (use_lake_model(i)>0) then
981 if (sncovr_ice(i) > f_zero) then
982 sfcemis_ice = emsref(7) * (f_one-sncovr_ice(i)) &
983 & + emsref(8) * sncovr_ice(i)
984 elseif (snodi(i) > f_zero) then
985 asnow = 0.02*snodi(i)
986 argh = min(0.50, max(.025,0.01*zorlf(i)))
987 hrgh = min(f_one,max(0.20,1.0577-1.1538e-3*hprif(i)))
988 fsno = asnow / (argh + asnow) * hrgh
989 sfcemis_ice = emsref(7)*(f_one-fsno) + emsref(8)*fsno
990 endif
991 semis_ice(i) = sfcemis_ice
992 else
993 sfcemis_ice = semis_ice(i) ! output from CICE or from RUC lsm (with snow effect)
994 endif
995 endif ! lsm check
996 endif ! icy
997
998 !-- land emissivity
999 !-- from Noah MP or RUC lsms
1000 sfcemis_land = semis_lnd(i) ! albedo with snow effect from LSM
1001
1002 !-- Composite emissivity from land, water and ice fractions.
1003 sfcemis(i) = fracl(i)*sfcemis_land + fraco(i)*emsref(1) &
1004 & + fraci(i)*sfcemis_ice
1005
1006 enddo ! i
1007
1008 endif ! end if_iemslw_block
1009
1010!chk print *,' In setemis, iemsflg, sfcemis =',iemsflg,sfcemis
1011
1012!
1013 return
1014!...................................
1015 end subroutine setemis
1016!-----------------------------------
1017
1018!.........................................!
1019 end module module_radiation_surface !
1021!=========================================!
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...
subroutine, public set_soilveg_ruc(me, isot, ivet, nlunit, errmsg, errflg)
This subroutine specifies vegetation and soil parameters for a given soil and land-use classification...
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, nf_albd, 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 the namelist options of soil/vegetation in RUC.
This module contains subroutine to specify vegetation and soil parameters for a given soild and land-...
This module contains routines used in the percentile matching algorithm for the albedo and vegetation...