albedo.F90 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277
  1. MODULE albedo
  2. !!======================================================================
  3. !! *** MODULE albedo ***
  4. !! Ocean forcing: bulk thermohaline forcing of the ocean (or ice)
  5. !!=====================================================================
  6. !! History : 8.0 ! 2001-04 (LIM 1.0)
  7. !! NEMO 1.0 ! 2003-07 (C. Ethe, G. Madec) Optimization (old name:shine)
  8. !! - ! 2004-11 (C. Talandier) add albedo_init
  9. !! - ! 2001-06 (M. Vancoppenolle) LIM 3.0
  10. !! - ! 2006-08 (G. Madec) cleaning for surface module
  11. !! 3.6 ! 2016-01 (C. Rousset) new parameterization for sea ice albedo
  12. !!----------------------------------------------------------------------
  13. !!----------------------------------------------------------------------
  14. !! albedo_ice : albedo for ice (clear and overcast skies)
  15. !! albedo_oce : albedo for ocean (clear and overcast skies)
  16. !! albedo_init : initialisation of albedo computation
  17. !!----------------------------------------------------------------------
  18. USE phycst ! physical constants
  19. USE in_out_manager ! I/O manager
  20. USE lib_mpp ! MPP library
  21. USE wrk_nemo ! work arrays
  22. USE lib_fortran ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined)
  23. IMPLICIT NONE
  24. PRIVATE
  25. PUBLIC albedo_ice ! routine called sbcice_lim.F90
  26. PUBLIC albedo_oce ! routine called by ???
  27. INTEGER :: albd_init = 0 !: control flag for initialization
  28. REAL(wp) :: rmue = 0.40 ! cosine of local solar altitude
  29. REAL(wp) :: ralb_oce = 0.066 ! ocean or lead albedo (Pegau and Paulson, Ann. Glac. 2001)
  30. REAL(wp) :: c1 = 0.05 ! snow thickness (only for nn_ice_alb=0)
  31. REAL(wp) :: c2 = 0.10 ! " "
  32. REAL(wp) :: rcloud = 0.06 ! cloud effect on albedo (only-for nn_ice_alb=0)
  33. ! !!* namelist namsbc_alb
  34. INTEGER :: nn_ice_alb
  35. REAL(wp) :: rn_alb_sdry, rn_alb_smlt, rn_alb_idry, rn_alb_imlt
  36. !!----------------------------------------------------------------------
  37. !! NEMO/OPA 3.3 , NEMO Consortium (2010)
  38. !! $Id: albedo.F90 4624 2014-04-28 12:09:03Z acc $
  39. !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
  40. !!----------------------------------------------------------------------
  41. CONTAINS
  42. SUBROUTINE albedo_ice( pt_ice, ph_ice, ph_snw, pa_ice_cs, pa_ice_os )
  43. !!----------------------------------------------------------------------
  44. !! *** ROUTINE albedo_ice ***
  45. !!
  46. !! ** Purpose : Computation of the albedo of the snow/ice system
  47. !!
  48. !! ** Method : Two schemes are available (from namelist parameter nn_ice_alb)
  49. !! 0: the scheme is that of Shine & Henderson-Sellers (JGR 1985) for clear-skies
  50. !! 1: the scheme is "home made" (for cloudy skies) and based on Brandt et al. (J. Climate 2005)
  51. !! and Grenfell & Perovich (JGR 2004)
  52. !! Description of scheme 1:
  53. !! 1) Albedo dependency on ice thickness follows the findings from Brandt et al (2005)
  54. !! which are an update of Allison et al. (JGR 1993) ; Brandt et al. 1999
  55. !! 0-5cm : linear function of ice thickness
  56. !! 5-150cm: log function of ice thickness
  57. !! > 150cm: constant
  58. !! 2) Albedo dependency on snow thickness follows the findings from Grenfell & Perovich (2004)
  59. !! i.e. it increases as -EXP(-snw_thick/0.02) during freezing and -EXP(-snw_thick/0.03) during melting
  60. !! 3) Albedo dependency on clouds is speculated from measurements of Grenfell and Perovich (2004)
  61. !! i.e. cloudy-clear albedo depend on cloudy albedo following a 2d order polynomial law
  62. !! 4) The needed 4 parameters are: dry and melting snow, freezing ice and bare puddled ice
  63. !!
  64. !! ** Note : The parameterization from Shine & Henderson-Sellers presents several misconstructions:
  65. !! 1) ice albedo when ice thick. tends to 0 is different than ocean albedo
  66. !! 2) for small ice thick. covered with some snow (<3cm?), albedo is larger
  67. !! under melting conditions than under freezing conditions
  68. !! 3) the evolution of ice albedo as a function of ice thickness shows
  69. !! 3 sharp inflexion points (at 5cm, 100cm and 150cm) that look highly unrealistic
  70. !!
  71. !! References : Shine & Henderson-Sellers 1985, JGR, 90(D1), 2243-2250.
  72. !! Brandt et al. 2005, J. Climate, vol 18
  73. !! Grenfell & Perovich 2004, JGR, vol 109
  74. !!----------------------------------------------------------------------
  75. REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: pt_ice ! ice surface temperature (Kelvin)
  76. REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: ph_ice ! sea-ice thickness
  77. REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: ph_snw ! snow thickness
  78. REAL(wp), INTENT( out), DIMENSION(:,:,:) :: pa_ice_cs ! albedo of ice under clear sky
  79. REAL(wp), INTENT( out), DIMENSION(:,:,:) :: pa_ice_os ! albedo of ice under overcast sky
  80. !!
  81. INTEGER :: ji, jj, jl ! dummy loop indices
  82. INTEGER :: ijpl ! number of ice categories (3rd dim of ice input arrays)
  83. REAL(wp) :: ralb_im, ralb_sf, ralb_sm, ralb_if
  84. REAL(wp) :: zswitch, z1_c1, z1_c2
  85. REAL(wp) :: zalb_sm, zalb_sf, zalb_st ! albedo of snow melting, freezing, total
  86. REAL(wp), POINTER, DIMENSION(:,:,:) :: zalb, zalb_it ! intermediate variable & albedo of ice (snow free)
  87. !!---------------------------------------------------------------------
  88. ijpl = SIZE( pt_ice, 3 ) ! number of ice categories
  89. CALL wrk_alloc( jpi,jpj,ijpl, zalb, zalb_it )
  90. IF( albd_init == 0 ) CALL albedo_init ! initialization
  91. ralb_sf = rn_alb_sdry ! dry snow
  92. ralb_sm = rn_alb_smlt ! melting snow
  93. ralb_if = rn_alb_idry ! bare frozen ice
  94. ralb_im = rn_alb_imlt ! bare puddled ice
  95. SELECT CASE ( nn_ice_alb )
  96. !------------------------------------------
  97. ! Shine and Henderson-Sellers (1985)
  98. !------------------------------------------
  99. CASE( 0 )
  100. ! Computation of ice albedo (free of snow)
  101. WHERE ( ph_snw == 0._wp .AND. pt_ice >= rt0_ice ) ; zalb(:,:,:) = ralb_im
  102. ELSE WHERE ; zalb(:,:,:) = ralb_if
  103. END WHERE
  104. WHERE ( 1.5 < ph_ice ) ; zalb_it = zalb
  105. ELSE WHERE( 1.0 < ph_ice .AND. ph_ice <= 1.5 ) ; zalb_it = 0.472 + 2.0 * ( zalb - 0.472 ) * ( ph_ice - 1.0 )
  106. ELSE WHERE( 0.05 < ph_ice .AND. ph_ice <= 1.0 ) ; zalb_it = 0.2467 + 0.7049 * ph_ice &
  107. & - 0.8608 * ph_ice * ph_ice &
  108. & + 0.3812 * ph_ice * ph_ice * ph_ice
  109. ELSE WHERE ; zalb_it = 0.1 + 3.6 * ph_ice
  110. END WHERE
  111. DO jl = 1, ijpl
  112. DO jj = 1, jpj
  113. DO ji = 1, jpi
  114. ! freezing snow
  115. ! no effect of underlying ice layer IF snow thickness > c1. Albedo does not depend on snow thick if > c2
  116. ! ! freezing snow
  117. zswitch = 1._wp - MAX( 0._wp , SIGN( 1._wp , - ( ph_snw(ji,jj,jl) - c1 ) ) )
  118. zalb_sf = ( 1._wp - zswitch ) * ( zalb_it(ji,jj,jl) &
  119. & + ph_snw(ji,jj,jl) * ( ralb_sf - zalb_it(ji,jj,jl) ) / c1 ) &
  120. & + zswitch * ralb_sf
  121. ! melting snow
  122. ! no effect of underlying ice layer. Albedo does not depend on snow thick IF > c2
  123. zswitch = MAX( 0._wp , SIGN( 1._wp , ph_snw(ji,jj,jl) - c2 ) )
  124. zalb_sm = ( 1._wp - zswitch ) * ( ralb_im + ph_snw(ji,jj,jl) * ( ralb_sm - ralb_im ) / c2 ) &
  125. & + zswitch * ralb_sm
  126. !
  127. ! snow albedo
  128. zswitch = MAX( 0._wp , SIGN( 1._wp , pt_ice(ji,jj,jl) - rt0_snow ) )
  129. zalb_st = zswitch * zalb_sm + ( 1._wp - zswitch ) * zalb_sf
  130. ! Ice/snow albedo
  131. zswitch = 1._wp - MAX( 0._wp , SIGN( 1._wp , - ph_snw(ji,jj,jl) ) )
  132. pa_ice_cs(ji,jj,jl) = zswitch * zalb_st + ( 1._wp - zswitch ) * zalb_it(ji,jj,jl)
  133. !
  134. END DO
  135. END DO
  136. END DO
  137. pa_ice_os(:,:,:) = pa_ice_cs(:,:,:) + rcloud ! Oberhuber correction for overcast sky
  138. !------------------------------------------
  139. ! New parameterization (2016)
  140. !------------------------------------------
  141. CASE( 1 )
  142. ! compilation of values from literature
  143. ! ralb_sf = 0.85 ! dry snow
  144. ! ralb_sm = 0.75 ! melting snow
  145. ! ralb_if = 0.60 ! bare frozen ice
  146. ! Perovich et al 2002 (Sheba) => the only dataset for which all types of ice/snow were retrieved
  147. ! ralb_sf = 0.85 ! dry snow
  148. ! ralb_sm = 0.72 ! melting snow
  149. ! ralb_if = 0.65 ! bare frozen ice
  150. ! Brandt et al 2005 (East Antarctica)
  151. ! ralb_sf = 0.87 ! dry snow
  152. ! ralb_sm = 0.82 ! melting snow
  153. ! ralb_if = 0.54 ! bare frozen ice
  154. !
  155. ! Computation of ice albedo (free of snow)
  156. z1_c1 = 1. / ( LOG(1.5) - LOG(0.05) )
  157. z1_c2 = 1. / 0.05
  158. WHERE ( ph_snw == 0._wp .AND. pt_ice >= rt0_ice ) ; zalb = ralb_im
  159. ELSE WHERE ; zalb = ralb_if
  160. END WHERE
  161. WHERE ( 1.5 < ph_ice ) ; zalb_it = zalb
  162. ELSE WHERE( 0.05 < ph_ice .AND. ph_ice <= 1.5 ) ; zalb_it = zalb + ( 0.18 - zalb ) * z1_c1 * &
  163. & ( LOG(1.5) - LOG(ph_ice) )
  164. ELSE WHERE ; zalb_it = ralb_oce + ( 0.18 - ralb_oce ) * z1_c2 * ph_ice
  165. END WHERE
  166. z1_c1 = 1. / 0.02
  167. z1_c2 = 1. / 0.03
  168. ! Computation of the snow/ice albedo
  169. DO jl = 1, ijpl
  170. DO jj = 1, jpj
  171. DO ji = 1, jpi
  172. zalb_sf = ralb_sf - ( ralb_sf - zalb_it(ji,jj,jl)) * EXP( - ph_snw(ji,jj,jl) * z1_c1 );
  173. zalb_sm = ralb_sm - ( ralb_sm - zalb_it(ji,jj,jl)) * EXP( - ph_snw(ji,jj,jl) * z1_c2 );
  174. ! snow albedo
  175. zswitch = MAX( 0._wp , SIGN( 1._wp , pt_ice(ji,jj,jl) - rt0_snow ) )
  176. zalb_st = zswitch * zalb_sm + ( 1._wp - zswitch ) * zalb_sf
  177. ! Ice/snow albedo
  178. zswitch = MAX( 0._wp , SIGN( 1._wp , - ph_snw(ji,jj,jl) ) )
  179. pa_ice_os(ji,jj,jl) = ( 1._wp - zswitch ) * zalb_st + zswitch * zalb_it(ji,jj,jl)
  180. END DO
  181. END DO
  182. END DO
  183. ! Effect of the clouds (2d order polynomial)
  184. pa_ice_cs = pa_ice_os - ( - 0.1010 * pa_ice_os * pa_ice_os + 0.1933 * pa_ice_os - 0.0148 );
  185. END SELECT
  186. CALL wrk_dealloc( jpi,jpj,ijpl, zalb, zalb_it )
  187. !
  188. END SUBROUTINE albedo_ice
  189. SUBROUTINE albedo_oce( pa_oce_os , pa_oce_cs )
  190. !!----------------------------------------------------------------------
  191. !! *** ROUTINE albedo_oce ***
  192. !!
  193. !! ** Purpose : Computation of the albedo of the ocean
  194. !!----------------------------------------------------------------------
  195. REAL(wp), DIMENSION(:,:), INTENT(out) :: pa_oce_os ! albedo of ocean under overcast sky
  196. REAL(wp), DIMENSION(:,:), INTENT(out) :: pa_oce_cs ! albedo of ocean under clear sky
  197. !!
  198. REAL(wp) :: zcoef
  199. !!----------------------------------------------------------------------
  200. !
  201. zcoef = 0.05 / ( 1.1 * rmue**1.4 + 0.15 ) ! Parameterization of Briegled and Ramanathan, 1982
  202. pa_oce_cs(:,:) = zcoef
  203. pa_oce_os(:,:) = 0.06 ! Parameterization of Kondratyev, 1969 and Payne, 1972
  204. !
  205. END SUBROUTINE albedo_oce
  206. SUBROUTINE albedo_init
  207. !!----------------------------------------------------------------------
  208. !! *** ROUTINE albedo_init ***
  209. !!
  210. !! ** Purpose : initializations for the albedo parameters
  211. !!
  212. !! ** Method : Read the namelist namsbc_alb
  213. !!----------------------------------------------------------------------
  214. INTEGER :: ios ! Local integer output status for namelist read
  215. NAMELIST/namsbc_alb/ nn_ice_alb, rn_alb_sdry, rn_alb_smlt, rn_alb_idry , rn_alb_imlt
  216. !!----------------------------------------------------------------------
  217. !
  218. albd_init = 1 ! indicate that the initialization has been done
  219. !
  220. REWIND( numnam_ref ) ! Namelist namsbc_alb in reference namelist : Albedo parameters
  221. READ ( numnam_ref, namsbc_alb, IOSTAT = ios, ERR = 901)
  222. 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_alb in reference namelist', lwp )
  223. REWIND( numnam_cfg ) ! Namelist namsbc_alb in configuration namelist : Albedo parameters
  224. READ ( numnam_cfg, namsbc_alb, IOSTAT = ios, ERR = 902 )
  225. 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_alb in configuration namelist', lwp )
  226. IF(lwm) WRITE ( numond, namsbc_alb )
  227. !
  228. IF(lwp) THEN ! Control print
  229. WRITE(numout,*)
  230. WRITE(numout,*) 'albedo : set albedo parameters'
  231. WRITE(numout,*) '~~~~~~~'
  232. WRITE(numout,*) ' Namelist namsbc_alb : albedo '
  233. WRITE(numout,*) ' choose the albedo parameterization nn_ice_alb = ', nn_ice_alb
  234. WRITE(numout,*) ' albedo of dry snow rn_alb_sdry = ', rn_alb_sdry
  235. WRITE(numout,*) ' albedo of melting snow rn_alb_smlt = ', rn_alb_smlt
  236. WRITE(numout,*) ' albedo of dry ice rn_alb_idry = ', rn_alb_idry
  237. WRITE(numout,*) ' albedo of bare puddled ice rn_alb_imlt = ', rn_alb_imlt
  238. ENDIF
  239. !
  240. END SUBROUTINE albedo_init
  241. !!======================================================================
  242. END MODULE albedo