Radiation Scheme in CCPP
module_radiation_surface Module Reference

This module sets up surface albedo for sw radiation and surface emissivity for lw radiation.

Functions/Subroutines

subroutine, public sfc_init
 This subroutine is the initialization program for surface radiation related quantities (albedo, emissivity, etc.) More...
 
subroutine, public setalb
 This subroutine computes four components of surface albedos (i.e., vis-nir, direct-diffused) according to control flag ialbflg.
1) climatological surface albedo scheme (briegleb 1992)
2) modis retrieval based scheme from boston univ. More...
 
subroutine, public setemis
 This subroutine computes surface emissivity for LW radiation. More...
 

Variables

character(40), parameter vtagsfc ='NCEP-Radiation_surface v5.1 Nov 2012 '
 
integer, parameter, public nf_albd = 4
 num of sfc albedo components More...
 
integer, parameter, public imxems = 360
 num of longitude points in global emis-type map More...
 
integer, parameter, public jmxems = 180
 num of latitude points in global emis-type map More...
 
real(kind=kind_phys), parameter f_zero = 0.0
 
real(kind=kind_phys), parameter f_one = 1.0
 
real(kind=kind_phys), parameter rad2dg = 180.0 / con_pi
 
integer, dimension(:,:), allocatable idxems
 global surface emissivity index array More...
 
integer iemslw = 0
 global surface emissivity contrl flag set up in 'sfc_init' More...
 

Function/Subroutine Documentation

subroutine, public module_radiation_surface::setalb ( )
Parameters
[in]slmskreal, (IMAX), sea(0),land(1),ice(2) mask on fcst model grid
[in]snowfreal, (IMAX), snow depth water equivalent in mm
[in]sncovrreal, (IMAX)
- ialgflg=0: not used
- ialgflg=1: snow cover over land in fraction
[in]snoalbreal, (IMAX)
- ialbflg=0: not used
- ialgflg=1: max snow albedo over land in fraction
[in]zorlfreal, (IMAX), surface roughness in cm
[in]coszfreal, (IMAX), cosin of solar zenith angle
[in]tsknfreal, (IMAX), ground surface temperature in K
[in]tairfreal, (IMAX), lowest model layer air temperature in K
[in]hprifreal, (IMAX), topographic sdv in m
— for ialbflg=0 climtological albedo scheme —
[in]alvsfreal, (IMAX), 60 degree vis albedo with strong cosz dependency
[in]alnsfreal, (IMAX), 60 degree nir albedo with strong cosz dependency
[in]alvwfreal, (IMAX), 60 degree vis albedo with weak cosz dependency
[in]alnwfreal, (IMAX), 60 degree nir albedo with weak cosz dependency
— for ialbflg=1 modis based land albedo scheme —
[in]alvsfreal, (IMAX), visible black sky albedo at zenith 60 degree
[in]alnsfreal, (IMAX), near-ir black sky albedo at zenith 60 degree
[in]alvwfreal, (IMAX), visible white sky albedo
[in]alnwfreal, (IMAX), near-ir white sky albedo
[in]facsfreal, (IMAX), fractional coverage with strong cosz dependency
[in]facwfreal, (IMAX), fractional coverage with weak cosz dependency
[in]ficereal, (IMAX), sea-ice fraction
[in]tisfcreal, (IMAX), sea-ice surface temperature
[in]IMAXinteger, array horizontal dimension
[out]sfcalbreal, (IMAX,NF_ALBD)
( :, 1) - near ir direct beam albedo
( :, 2) - near ir diffused albedo
( :, 3) - uv+vis direct beam albedo
( :, 4) - uv+vis diffused albedo

General Algorithm

  1. If use climatological albedo scheme:
    • Modified snow albedo scheme - units convert to m (originally snowf in mm; zorlf in cm)
    • Compute diffused sea surface albedo
    • Compute diffused snow albedo
    • Compute direct snow albedo
    • Compute direct sea surface albedo
  2. If use modis based albedo for land area:
    • Compute snow cover input directly for land model, no conversion needed.
    • Compute diffused sea surface albedo
    • Compute diffused snow albedo, land area use input max snow albedo
    • Compute direct snow albedo
    • Compute direct sea surface albedo, use fanglin's zenith angle treatment

Definition at line 303 of file radiation_surface.f.

References physcons::con_t0c, physcons::con_tice, physcons::con_ttp, f_one, f_zero, and physparam::ialbflg.

Referenced by module_radiation_driver::grrad().

303 
304 ! --- inputs:
305  & ( slmsk,snowf,sncovr,snoalb,zorlf,coszf,tsknf,tairf,hprif, &
306  & alvsf,alnsf,alvwf,alnwf,facsf,facwf,fice,tisfc, &
307  & imax, &
308 ! --- outputs:
309  & sfcalb &
310  & )
311 
312 ! =================================================================== !
313 ! !
314 ! this program computes four components of surface albedos (i.e. !
315 ! vis-nir, direct-diffused) according to controflag ialbflg. !
316 ! 1) climatological surface albedo scheme (briegleb 1992) !
317 ! 2) modis retrieval based scheme from boston univ. !
318 ! !
319 ! !
320 ! usage: call setalb !
321 ! !
322 ! subprograms called: none !
323 ! !
324 ! ==================== defination of variables ==================== !
325 ! !
326 ! inputs: !
327 ! slmsk (IMAX) - sea(0),land(1),ice(2) mask on fcst model grid !
328 ! snowf (IMAX) - snow depth water equivalent in mm !
329 ! sncovr(IMAX) - ialgflg=0: not used !
330 ! ialgflg=1: snow cover over land in fraction !
331 ! snoalb(IMAX) - ialbflg=0: not used !
332 ! ialgflg=1: max snow albedo over land in fraction !
333 ! zorlf (IMAX) - surface roughness in cm !
334 ! coszf (IMAX) - cosin of solar zenith angle !
335 ! tsknf (IMAX) - ground surface temperature in k !
336 ! tairf (IMAX) - lowest model layer air temperature in k !
337 ! hprif (IMAX) - topographic sdv in m !
338 ! --- for ialbflg=0 climtological albedo scheme --- !
339 ! alvsf (IMAX) - 60 degree vis albedo with strong cosz dependency !
340 ! alnsf (IMAX) - 60 degree nir albedo with strong cosz dependency !
341 ! alvwf (IMAX) - 60 degree vis albedo with weak cosz dependency !
342 ! alnwf (IMAX) - 60 degree nir albedo with weak cosz dependency !
343 ! --- for ialbflg=1 modis based land albedo scheme --- !
344 ! alvsf (IMAX) - visible black sky albedo at zenith 60 degree !
345 ! alnsf (IMAX) - near-ir black sky albedo at zenith 60 degree !
346 ! alvwf (IMAX) - visible white sky albedo !
347 ! alnwf (IMAX) - near-ir white sky albedo !
348 ! !
349 ! facsf (IMAX) - fractional coverage with strong cosz dependency !
350 ! facwf (IMAX) - fractional coverage with weak cosz dependency !
351 ! fice (IMAX) - sea-ice fraction !
352 ! tisfc (IMAX) - sea-ice surface temperature !
353 ! IMAX - array horizontal dimension !
354 ! !
355 ! outputs: !
356 ! sfcalb(IMAX,NF_ALBD) !
357 ! ( :, 1) - near ir direct beam albedo !
358 ! ( :, 2) - near ir diffused albedo !
359 ! ( :, 3) - uv+vis direct beam albedo !
360 ! ( :, 4) - uv+vis diffused albedo !
361 ! !
362 ! module internal control variables: !
363 ! ialbflg - =0 use the default climatology surface albedo !
364 ! =1 use modis retrieved albedo and input snow cover!
365 ! for land areas !
366 ! !
367 ! ==================== end of description ===================== !
368 !
369  implicit none
370 
371 ! --- inputs
372  integer, intent(in) :: imax
373 
374  real (kind=kind_phys), dimension(:), intent(in) :: &
375  & slmsk, snowf, zorlf, coszf, tsknf, tairf, hprif, &
376  & alvsf, alnsf, alvwf, alnwf, facsf, facwf, fice, tisfc, &
377  & sncovr, snoalb
378 
379 ! --- outputs
380  real (kind=kind_phys), dimension(IMAX,NF_ALBD), intent(out) :: &
381  & sfcalb
382 ! real (kind=kind_phys), dimension(:,:), intent(out) :: sfcalb
383 
384 ! --- locals:
385  real (kind=kind_phys) :: asnvb, asnnb, asnvd, asnnd, asevb &
386  &, asenb, asevd, asend, fsno, fsea, rfcs, rfcw, flnd &
387  &, asnow, argh, hrgh, fsno0, fsno1, flnd0, fsea0, csnow &
388  &, a1, a2, b1, b2, b3, ab1bm, ab2bm
389 
390  real (kind=kind_phys) ffw, dtgd
391 
392  integer :: i, k
393 
394 !
395 !===> ... begin here
396 !
398  if ( ialbflg == 0 ) then ! use climatological albedo scheme
399 
400  do i = 1, imax
401 
403 ! --- modified snow albedo scheme - units convert to m
404 ! (originally snowf in mm; zorlf in cm)
405 
406  asnow = 0.02*snowf(i)
407  argh = min(0.50, max(.025, 0.01*zorlf(i)))
408  hrgh = min(f_one, max(0.20, 1.0577-1.1538e-3*hprif(i) ) )
409  fsno0 = asnow / (argh + asnow) * hrgh
410  if (nint(slmsk(i))==0 .and. tsknf(i)>con_tice) fsno0 = f_zero
411  fsno1 = f_one - fsno0
412  flnd0 = min(f_one, facsf(i)+facwf(i))
413  fsea0 = max(f_zero, f_one-flnd0)
414  fsno = fsno0
415  fsea = fsea0 * fsno1
416  flnd = flnd0 * fsno1
417 
419 ! --- diffused sea surface albedo
420 
421  if (tsknf(i) >= 271.5) then
422  asevd = 0.06
423  asend = 0.06
424  elseif (tsknf(i) < 271.1) then
425  asevd = 0.70
426  asend = 0.65
427  else
428  a1 = (tsknf(i) - 271.1)**2
429  asevd = 0.7 - 4.0*a1
430  asend = 0.65 - 3.6875*a1
431  endif
432 
434 ! --- diffused snow albedo
435 
436  if (nint(slmsk(i)) == 2) then
437  ffw = f_one - fice(i)
438  if (ffw < f_one) then
439  dtgd = max(f_zero, min(5.0, (con_ttp-tisfc(i)) ))
440  b1 = 0.03 * dtgd
441  else
442  b1 = f_zero
443  endif
444 
445  b3 = 0.06 * ffw
446  asnvd = (0.70 + b1) * fice(i) + b3
447  asnnd = (0.60 + b1) * fice(i) + b3
448  asevd = 0.70 * fice(i) + b3
449  asend = 0.60 * fice(i) + b3
450  else
451  asnvd = 0.90
452  asnnd = 0.75
453  endif
454 
456 ! --- direct snow albedo
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 ! --- direct sea surface albedo
469 
470  if (coszf(i) > 0.0001) then
471 ! rfcs = 1.4 / (f_one + 0.8*coszf(i))
472 ! rfcw = 1.3 / (f_one + 0.6*coszf(i))
473  rfcs = 2.14 / (f_one + 1.48*coszf(i))
474  rfcw = rfcs
475 
476  if (tsknf(i) >= con_t0c) then
477  asevb = max(asevd, 0.026/(coszf(i)**1.7+0.065) &
478  & + 0.15 * (coszf(i)-0.1) * (coszf(i)-0.5) &
479  & * (coszf(i)-f_one))
480  asenb = asevb
481  else
482  asevb = asevd
483  asenb = asend
484  endif
485  else
486  rfcs = f_one
487  rfcw = f_one
488  asevb = asevd
489  asenb = asend
490  endif
491 
492  a1 = alvsf(i) * facsf(i)
493  b1 = alvwf(i) * facwf(i)
494  a2 = alnsf(i) * facsf(i)
495  b2 = alnwf(i) * facwf(i)
496  ab1bm = a1*rfcs + b1*rfcw
497  ab2bm = a2*rfcs + b2*rfcw
498  sfcalb(i,1) = min(0.99, ab2bm) *flnd + asenb*fsea + asnnb*fsno
499  sfcalb(i,2) = (a2 + b2) * 0.96 *flnd + asend*fsea + asnnd*fsno
500  sfcalb(i,3) = min(0.99, ab1bm) *flnd + asevb*fsea + asnvb*fsno
501  sfcalb(i,4) = (a1 + b1) * 0.96 *flnd + asevd*fsea + asnvd*fsno
502 
503  enddo ! end_do_i_loop
504 
506  else ! use modis based albedo for land area
507 
508  do i = 1, imax
509 
511 ! --- snow cover input directly from land model, no conversion needed
512 
513  fsno0 = sncovr(i)
514 
515  if (nint(slmsk(i))==0 .and. tsknf(i)>con_tice) fsno0 = f_zero
516 
517  if (nint(slmsk(i)) == 2) then
518  asnow = 0.02*snowf(i)
519  argh = min(0.50, max(.025, 0.01*zorlf(i)))
520  hrgh = min(f_one, max(0.20, 1.0577-1.1538e-3*hprif(i) ) )
521  fsno0 = asnow / (argh + asnow) * hrgh
522  endif
523 
524  fsno1 = f_one - fsno0
525  flnd0 = min(f_one, facsf(i)+facwf(i))
526  fsea0 = max(f_zero, f_one-flnd0)
527  fsno = fsno0
528  fsea = fsea0 * fsno1
529  flnd = flnd0 * fsno1
530 
532 ! --- diffused sea surface albedo
533 
534  if (tsknf(i) >= 271.5) then
535  asevd = 0.06
536  asend = 0.06
537  elseif (tsknf(i) < 271.1) then
538  asevd = 0.70
539  asend = 0.65
540  else
541  a1 = (tsknf(i) - 271.1)**2
542  asevd = 0.7 - 4.0*a1
543  asend = 0.65 - 3.6875*a1
544  endif
545 
547 ! --- diffused snow albedo, land area use input max snow albedo
548 
549  if (nint(slmsk(i)) == 2) then
550  ffw = f_one - fice(i)
551  if (ffw < f_one) then
552  dtgd = max(f_zero, min(5.0, (con_ttp-tisfc(i)) ))
553  b1 = 0.03 * dtgd
554  else
555  b1 = f_zero
556  endif
557 
558  b3 = 0.06 * ffw
559  asnvd = (0.70 + b1) * fice(i) + b3
560  asnnd = (0.60 + b1) * fice(i) + b3
561  asevd = 0.70 * fice(i) + b3
562  asend = 0.60 * fice(i) + b3
563  else
564  asnvd = snoalb(i)
565  asnnd = snoalb(i)
566  endif
567 
569 ! --- direct snow albedo
570 
571  if (nint(slmsk(i)) == 2) then
572  if (coszf(i) < 0.5) then
573  csnow = 0.5 * (3.0 / (f_one+4.0*coszf(i)) - f_one)
574  asnvb = min( 0.98, asnvd+(f_one-asnvd)*csnow )
575  asnnb = min( 0.98, asnnd+(f_one-asnnd)*csnow )
576  else
577  asnvb = asnvd
578  asnnb = asnnd
579  endif
580  else
581  asnvb = snoalb(i)
582  asnnb = snoalb(i)
583  endif
584 
586 ! --- direct sea surface albedo, use fanglin's zenith angle treatment
587 
588  if (coszf(i) > 0.0001) then
589 
590 ! rfcs = 1.89 - 3.34*coszf(i) + 4.13*coszf(i)*coszf(i) &
591 ! & - 2.02*coszf(i)*coszf(i)*coszf(i)
592  rfcs = 1.775/(1.0+1.55*coszf(i))
593 
594  if (tsknf(i) >= con_t0c) then
595  asevb = max(asevd, 0.026/(coszf(i)**1.7+0.065) &
596  & + 0.15 * (coszf(i)-0.1) * (coszf(i)-0.5) &
597  & * (coszf(i)-f_one))
598  asenb = asevb
599  else
600  asevb = asevd
601  asenb = asend
602  endif
603  else
604  rfcs = f_one
605  asevb = asevd
606  asenb = asend
607  endif
608 
609  ab1bm = min(0.99, alnsf(i)*rfcs)
610  ab2bm = min(0.99, alvsf(i)*rfcs)
611  sfcalb(i,1) = ab1bm *flnd + asenb*fsea + asnnb*fsno
612  sfcalb(i,2) = alnwf(i) *flnd + asend*fsea + asnnd*fsno
613  sfcalb(i,3) = ab2bm *flnd + asevb*fsea + asnvb*fsno
614  sfcalb(i,4) = alvwf(i) *flnd + asevd*fsea + asnvd*fsno
615 
616  enddo ! end_do_i_loop
617 
618  endif ! end if_ialbflg
619 !
620  return
621 !...................................

Here is the caller graph for this function:

subroutine, public module_radiation_surface::setemis ( )
Parameters
[in]xlonreal, (IMAX), longitude in radiance, ok for both 0->2pi or -pi -> +pi ranges
[in]xlatreal, (IMAX), latitude in radiance, default to pi/2 -> -pi/2 range, otherwise see in-line comment
[in]slmskreal, (IMAX), sea(0),land(1),ice(2) mask on fcst model grid
[in]snowfreal, (IMAX), snow depth water equivalent in mm
[in]sncovrreal, (IMAX), ialbflg=1: snow cover over land in fraction
[in]zorlfreal, (IMAX), surface roughness in cm
[in]tsknfreal, (IMAX), ground surface temperature in K
[in]tairfreal, (IMAX), lowest model layer air temperature in K
[in]hprifreal, (IMAX), topographic sdv in m
[in]IMAXinteger, array horizontal dimension
[out]sfcemisreal, (IMAX), surface emissivity

General Algorithm

  1. Map input data onto model grid
  2. Map grid in longitude direction
  3. Map grid in latitude direction
  4. Check for snow covered area

Definition at line 643 of file radiation_surface.f.

References f_one, f_zero, physparam::ialbflg, idxems, iemslw, imxems, jmxems, and rad2dg.

Referenced by module_radiation_driver::grrad().

643 
644 ! --- inputs:
645  & ( xlon,xlat,slmsk,snowf,sncovr,zorlf,tsknf,tairf,hprif, &
646  & imax, &
647 ! --- outputs:
648  & sfcemis &
649  & )
650 
651 ! =================================================================== !
652 ! !
653 ! this program computes surface emissivity for lw radiation. !
654 ! !
655 ! usage: call setemis !
656 ! !
657 ! subprograms called: none !
658 ! !
659 ! ==================== defination of variables ==================== !
660 ! !
661 ! inputs: !
662 ! xlon (IMAX) - longitude in radiance, ok for both 0->2pi or !
663 ! -pi -> +pi ranges !
664 ! xlat (IMAX) - latitude in radiance, default to pi/2 -> -pi/2 !
665 ! range, otherwise see in-line comment !
666 ! slmsk (IMAX) - sea(0),land(1),ice(2) mask on fcst model grid !
667 ! snowf (IMAX) - snow depth water equivalent in mm !
668 ! sncovr(IMAX) - ialbflg=1: snow cover over land in fraction !
669 ! zorlf (IMAX) - surface roughness in cm !
670 ! tsknf (IMAX) - ground surface temperature in k !
671 ! tairf (IMAX) - lowest model layer air temperature in k !
672 ! hprif (IMAX) - topographic sdv in m !
673 ! IMAX - array horizontal dimension !
674 ! !
675 ! outputs: !
676 ! sfcemis(IMAX) - surface emissivity !
677 ! !
678 ! ------------------------------------------------------------------- !
679 ! !
680 ! surface type definations: !
681 ! 1. open water 2. grass/wood/shrub land !
682 ! 3. tundra/bare soil 4. sandy desert !
683 ! 5. rocky desert 6. forest !
684 ! 7. ice 8. snow !
685 ! !
686 ! input index data lon from 0 towards east, lat from n to s !
687 ! !
688 ! ==================== end of description ===================== !
689 !
690  implicit none
691 
692 ! --- inputs
693  integer, intent(in) :: imax
694 
695  real (kind=kind_phys), dimension(:), intent(in) :: &
696  & xlon,xlat, slmsk, snowf,sncovr, zorlf, tsknf, tairf, hprif
697 
698 ! --- outputs
699  real (kind=kind_phys), dimension(:), intent(out) :: sfcemis
700 
701 ! --- locals:
702  integer :: i, i1, i2, j1, j2, idx
703 
704  real (kind=kind_phys) :: dltg, hdlt, tmp1, tmp2, &
705  & asnow, argh, hrgh, fsno, fsno0, fsno1
706 
707 ! --- reference emiss value for diff surface emiss index
708 ! 1-open water, 2-grass/shrub land, 3-bare soil, tundra,
709 ! 4-sandy desert, 5-rocky desert, 6-forest, 7-ice, 8-snow
710 
711  real (kind=kind_phys) :: emsref(8)
712  data emsref / 0.97, 0.95, 0.94, 0.90, 0.93, 0.96, 0.96, 0.99 /
713 
714 !
715 !===> ... begin here
716 !
717 
718  if ( iemslw == 0 ) then ! sfc emiss default to 1.0
719 
720  sfcemis(:) = f_one
721  return
722 
723  else ! emiss set by sfc type and condition
724 
725  dltg = 360.0 / float(imxems)
726  hdlt = 0.5 * dltg
727 
729 ! --- ... mapping input data onto model grid
730 ! note: this is a simple mapping method, an upgrade is needed if
731 ! the model grid is much corcer than the 1-deg data resolution
732 
733  lab_do_imax : do i = 1, imax
734 
735  if ( nint(slmsk(i)) == 0 ) then ! sea point
736 
737  sfcemis(i) = emsref(1)
738 
739  else if ( nint(slmsk(i)) == 2 ) then ! sea-ice
740 
741  sfcemis(i) = emsref(7)
742 
743  else ! land
745 ! --- map grid in longitude direction
746  i2 = 1
747  j2 = 1
748 
749  tmp1 = xlon(i) * rad2dg
750  if (tmp1 < f_zero) tmp1 = tmp1 + 360.0
751 
752  lab_do_imxems : do i1 = 1, imxems
753  tmp2 = dltg * (i1 - 1) + hdlt
754 
755  if (abs(tmp1-tmp2) <= hdlt) then
756  i2 = i1
757  exit lab_do_imxems
758  endif
759  enddo lab_do_imxems
761 ! --- map grid in latitude direction
762  tmp1 = xlat(i) * rad2dg ! if xlat in pi/2 -> -pi/2 range
763 ! tmp1 = 90.0 - xlat(i)*rad2dg ! if xlat in 0 -> pi range
764 
765  lab_do_jmxems : do j1 = 1, jmxems
766  tmp2 = 90.0 - dltg * (j1 - 1)
767 
768  if (abs(tmp1-tmp2) <= hdlt) then
769  j2 = j1
770  exit lab_do_jmxems
771  endif
772  enddo lab_do_jmxems
773 
774 
775  idx = max( 2, idxems(i2,j2) )
776  if ( idx >= 7 ) idx = 2
777  sfcemis(i) = emsref(idx)
778 
779  endif ! end if_slmsk_block
780 
782 ! --- check for snow covered area
783 
784  if ( ialbflg==1 .and. nint(slmsk(i))==1 ) then ! input land area snow cover
785 
786  fsno0 = sncovr(i)
787  fsno1 = f_one - fsno0
788  sfcemis(i) = sfcemis(i)*fsno1 + emsref(8)*fsno0
789 
790  else ! compute snow cover from snow depth
791  if ( snowf(i) > f_zero ) then
792  asnow = 0.02*snowf(i)
793  argh = min(0.50, max(.025, 0.01*zorlf(i)))
794  hrgh = min(f_one, max(0.20, 1.0577-1.1538e-3*hprif(i) ) )
795  fsno0 = asnow / (argh + asnow) * hrgh
796  if (nint(slmsk(i)) == 0 .and. tsknf(i) > 271.2) &
797  & fsno0=f_zero
798  fsno1 = f_one - fsno0
799  sfcemis(i) = sfcemis(i)*fsno1 + emsref(8)*fsno0
800  endif
801 
802  endif ! end if_ialbflg
803 
804  enddo lab_do_imax
805 
806  endif ! end if_iemslw_block
807 
808 !chk print *,' In setemis, iemsflg, sfcemis =',iemsflg,sfcemis
809 
810 !
811  return
812 !...................................

Here is the caller graph for this function:

subroutine, public module_radiation_surface::sfc_init ( )
Parameters
[in]meinteger, print control flag
[out]NONE

External Module Variables


physparam::ialbflg - control flag for surface albedo schemes.
=0: climatology, based on surface veg types; =1:
physparam::iemsflg - control flag for sfc emissivity schemes (ab:2-dig)
a:=0 set sfc air/ground t same for lw radiation; =1 set sfc air/ground t diff for lw radiation
b:=0 use fixed sfc emissivity=1.0 (black-body); =1 use varying climtology sfc emiss (veg based)

Definition at line 129 of file radiation_surface.f.

References physparam::ialbflg, idxems, physparam::iemsflg, iemslw, imxems, jmxems, physparam::semis_file, and vtagsfc.

Referenced by module_radiation_driver::radinit().

129 
130 ! --- inputs:
131  & ( me )
132 ! --- outputs: ( none )
133 
134 ! =================================================================== !
135 ! !
136 ! this program is the initialization program for surface radiation !
137 ! related quantities (albedo, emissivity, etc.) !
138 ! !
139 ! usage: call sfc_init !
140 ! !
141 ! subprograms called: none !
142 ! !
143 ! ==================== defination of variables ==================== !
144 ! !
145 ! inputs: !
146 ! me - print control flag !
147 ! !
148 ! outputs: (none) to module variables only !
149 ! !
150 ! external module variables: !
151 ! ialbflg - control flag for surface albedo schemes !
152 ! =0: climatology, based on surface veg types !
153 ! =1: !
154 ! iemsflg - control flag for sfc emissivity schemes (ab:2-dig)!
155 ! a:=0 set sfc air/ground t same for lw radiation !
156 ! =1 set sfc air/ground t diff for lw radiation !
157 ! b:=0 use fixed sfc emissivity=1.0 (black-body) !
158 ! =1 use varying climtology sfc emiss (veg based) !
159 ! !
160 ! ==================== end of description ===================== !
161 !
162  implicit none
163 
164 ! --- inputs:
165  integer, intent(in) :: me
166 
167 ! --- outputs: ( none )
168 
169 ! --- locals:
170  integer :: i, k
171 ! integer :: ia, ja
172  logical :: file_exist
173  character :: cline*80
174 !
175 !===> ... begin here
176 !
177  if ( me == 0 ) print *, vtagsfc ! print out version tag
178 
179 ! --- ... initialization of surface albedo section
180 
181  if ( ialbflg == 0 ) then
182 
183  if ( me == 0 ) then
184  print *,' - Using climatology surface albedo scheme for sw'
185  endif
186 
187  else if ( ialbflg == 1 ) then
188 
189  if ( me == 0 ) then
190  print *,' - Using MODIS based land surface albedo for sw'
191  endif
192 
193  else
194  print *,' !! ERROR in Albedo Scheme Setting, IALB=',ialbflg
195  stop
196  endif ! end if_ialbflg_block
197 
198 ! --- ... initialization of surface emissivity section
199 
200  iemslw = mod(iemsflg, 10) ! emissivity control
201  if ( iemslw == 0 ) then ! fixed sfc emis at 1.0
202 
203  if ( me == 0 ) then
204  print *,' - Using Fixed Surface Emissivity = 1.0 for lw'
205  endif
206 
207  elseif ( iemslw == 1 ) then ! input sfc emiss type map
208 
209 ! --- allocate data space
210  if ( .not. allocated(idxems) ) then
211  allocate ( idxems(imxems,jmxems) )
212  endif
213 
214 ! --- check to see if requested emissivity data file existed
215 
216  inquire (file=semis_file, exist=file_exist)
217 
218  if ( .not. file_exist ) then
219  if ( me == 0 ) then
220  print *,' - Using Varying Surface Emissivity for lw'
221  print *,' Requested data file "',semis_file,'" not found!'
222  print *,' Change to fixed surface emissivity = 1.0 !'
223  endif
224 
225  iemslw = 0
226  else
227  close(niradsf)
228  open (niradsf,file=semis_file,form='formatted',status='old')
229  rewind niradsf
230 
231  read (niradsf,12) cline
232  12 format(a80)
233 
234  read (niradsf,14) idxems
235  14 format(80i1)
236 
237  if ( me == 0 ) then
238  print *,' - Using Varying Surface Emissivity for lw'
239  print *,' Opened data file: ',semis_file
240  print *, cline
241 !check print *,' CHECK: Sample emissivity index data'
242 ! ia = IMXEMS / 5
243 ! ja = JMXEMS / 5
244 ! print *, idxems(1:IMXEMS:ia,1:JMXEMS:ja)
245  endif
246 
247  close(niradsf)
248  endif ! end if_file_exist_block
249 
250  else
251  print *,' !! ERROR in Emissivity Scheme Setting, IEMS=',iemsflg
252  stop
253  endif ! end if_iemslw_block
254 
255 !
256  return
257 !...................................

Here is the caller graph for this function:

Variable Documentation

real (kind=kind_phys), parameter module_radiation_surface::f_one = 1.0
private

Definition at line 101 of file radiation_surface.f.

Referenced by setalb(), and setemis().

101  real (kind=kind_phys), parameter :: f_one = 1.0
real (kind=kind_phys), parameter module_radiation_surface::f_zero = 0.0
private

Definition at line 100 of file radiation_surface.f.

Referenced by setalb(), and setemis().

100  real (kind=kind_phys), parameter :: f_zero = 0.0
integer, dimension(:,:), allocatable module_radiation_surface::idxems
private

Definition at line 106 of file radiation_surface.f.

Referenced by setemis(), and sfc_init().

106  integer, allocatable :: idxems(:,:)
integer module_radiation_surface::iemslw = 0
private

Definition at line 108 of file radiation_surface.f.

Referenced by setemis(), and sfc_init().

108  integer :: iemslw = 0
integer, parameter, public module_radiation_surface::imxems = 360

Definition at line 96 of file radiation_surface.f.

Referenced by setemis(), and sfc_init().

96  integer, parameter, public :: imxems = 360 ! num of lon-pts in glb emis-type map
integer, parameter, public module_radiation_surface::jmxems = 180

Definition at line 98 of file radiation_surface.f.

Referenced by setemis(), and sfc_init().

98  integer, parameter, public :: jmxems = 180 ! num of lat-pts in glb emis-type map
integer, parameter, public module_radiation_surface::nf_albd = 4

Definition at line 94 of file radiation_surface.f.

94  integer, parameter, public :: nf_albd = 4 ! num of sfc albedo components
real (kind=kind_phys), parameter module_radiation_surface::rad2dg = 180.0 / con_pi
private

Definition at line 102 of file radiation_surface.f.

Referenced by setemis().

102  real (kind=kind_phys), parameter :: rad2dg= 180.0 / con_pi
character(40), parameter module_radiation_surface::vtagsfc ='NCEP-Radiation_surface v5.1 Nov 2012 '
private

Definition at line 88 of file radiation_surface.f.

Referenced by sfc_init().

88  character(40), parameter :: &
89  & VTAGSFC='NCEP-Radiation_surface v5.1 Nov 2012 '