1 | MODULE icealb |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE icealb *** |
---|
4 | !! Atmospheric forcing: Albedo over sea ice |
---|
5 | !!===================================================================== |
---|
6 | !! History : 4.0 ! 2017-07 (C. Rousset) Split ice and ocean albedos |
---|
7 | !! 4.0 ! 2018 (many people) SI3 [aka Sea Ice cube] |
---|
8 | !!---------------------------------------------------------------------- |
---|
9 | #if defined key_si3 |
---|
10 | !!---------------------------------------------------------------------- |
---|
11 | !! 'key_si3' SI3 sea-ice model |
---|
12 | !!---------------------------------------------------------------------- |
---|
13 | !! ice_alb : albedo for ice (clear and overcast skies) |
---|
14 | !! ice_alb_init : initialisation of albedo computation |
---|
15 | !!---------------------------------------------------------------------- |
---|
16 | USE ice, ONLY: jpl ! sea-ice: number of categories |
---|
17 | USE phycst ! physical constants |
---|
18 | USE dom_oce ! domain: ocean |
---|
19 | ! |
---|
20 | USE in_out_manager ! I/O manager |
---|
21 | USE lib_mpp ! MPP library |
---|
22 | USE lib_fortran ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) |
---|
23 | USE timing ! Timing |
---|
24 | |
---|
25 | IMPLICIT NONE |
---|
26 | PRIVATE |
---|
27 | |
---|
28 | PUBLIC ice_alb_init ! called in icestp |
---|
29 | PUBLIC ice_alb ! called in icesbc.F90 and iceupdate.F90 |
---|
30 | |
---|
31 | REAL(wp), PUBLIC, PARAMETER :: rn_alb_oce = 0.066 !: ocean or lead albedo (Pegau and Paulson, Ann. Glac. 2001) |
---|
32 | ! |
---|
33 | ! !!* albedo namelist (namalb) |
---|
34 | REAL(wp) :: rn_alb_sdry ! dry snow albedo |
---|
35 | REAL(wp) :: rn_alb_smlt ! melting snow albedo |
---|
36 | REAL(wp) :: rn_alb_idry ! dry ice albedo |
---|
37 | REAL(wp) :: rn_alb_imlt ! bare puddled ice albedo |
---|
38 | REAL(wp) :: rn_alb_dpnd ! ponded ice albedo |
---|
39 | |
---|
40 | !!---------------------------------------------------------------------- |
---|
41 | !! NEMO/ICE 4.0 , NEMO Consortium (2018) |
---|
42 | !! $Id$ |
---|
43 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
44 | !!---------------------------------------------------------------------- |
---|
45 | CONTAINS |
---|
46 | |
---|
47 | SUBROUTINE ice_alb( pt_su, ph_ice, ph_snw, ld_pnd_alb, pafrac_pnd, ph_pnd, plh_pnd, palb_cs, palb_os ) |
---|
48 | !!---------------------------------------------------------------------- |
---|
49 | !! *** ROUTINE ice_alb *** |
---|
50 | !! |
---|
51 | !! ** Purpose : Computation of the albedo of the snow/ice system |
---|
52 | !! |
---|
53 | !! ** Method : The scheme is "home made" (for cloudy skies) and based on Brandt et al. (J. Climate 2005) |
---|
54 | !! and Grenfell & Perovich (JGR 2004) |
---|
55 | !! 1) Albedo dependency on ice thickness follows the findings from Brandt et al (2005) |
---|
56 | !! which are an update of Allison et al. (JGR 1993) ; Brandt et al. 1999 |
---|
57 | !! 0-5cm : linear function of ice thickness |
---|
58 | !! 5-150cm: log function of ice thickness |
---|
59 | !! > 150cm: constant |
---|
60 | !! 2) Albedo dependency on snow thickness follows the findings from Grenfell & Perovich (2004) |
---|
61 | !! i.e. it increases as -EXP(-snw_thick/0.02) during freezing and -EXP(-snw_thick/0.03) during melting |
---|
62 | !! 3) Albedo dependency on clouds is speculated from measurements of Grenfell and Perovich (2004) |
---|
63 | !! i.e. cloudy-clear albedo depend on cloudy albedo following a 2d order polynomial law |
---|
64 | !! 4) The needed 4 parameters are: dry and melting snow, freezing ice and bare puddled ice |
---|
65 | !! |
---|
66 | !! compilation of values from literature (reference overcast sky values) |
---|
67 | !! rn_alb_sdry = 0.85 ! dry snow |
---|
68 | !! rn_alb_smlt = 0.75 ! melting snow |
---|
69 | !! rn_alb_idry = 0.60 ! bare frozen ice |
---|
70 | !! rn_alb_imlt = 0.50 ! bare puddled ice albedo |
---|
71 | !! rn_alb_dpnd = 0.36 ! ponded-ice overcast albedo (Lecomte et al, 2015) |
---|
72 | !! ! early melt pnds 0.27, late melt ponds 0.14 Grenfell & Perovich |
---|
73 | !! Perovich et al 2002 (Sheba) => the only dataset for which all types of ice/snow were retrieved |
---|
74 | !! rn_alb_sdry = 0.85 ! dry snow |
---|
75 | !! rn_alb_smlt = 0.72 ! melting snow |
---|
76 | !! rn_alb_idry = 0.65 ! bare frozen ice |
---|
77 | !! Brandt et al 2005 (East Antarctica) |
---|
78 | !! rn_alb_sdry = 0.87 ! dry snow |
---|
79 | !! rn_alb_smlt = 0.82 ! melting snow |
---|
80 | !! rn_alb_idry = 0.54 ! bare frozen ice |
---|
81 | !! |
---|
82 | !! ** Note : The old parameterization from Shine & Henderson-Sellers (not here anymore) presented several misconstructions |
---|
83 | !! 1) ice albedo when ice thick. tends to 0 is different than ocean albedo |
---|
84 | !! 2) for small ice thick. covered with some snow (<3cm?), albedo is larger |
---|
85 | !! under melting conditions than under freezing conditions |
---|
86 | !! 3) the evolution of ice albedo as a function of ice thickness shows |
---|
87 | !! 3 sharp inflexion points (at 5cm, 100cm and 150cm) that look highly unrealistic |
---|
88 | !! |
---|
89 | !! References : Shine & Henderson-Sellers 1985, JGR, 90(D1), 2243-2250. |
---|
90 | !! Brandt et al. 2005, J. Climate, vol 18 |
---|
91 | !! Grenfell & Perovich 2004, JGR, vol 109 |
---|
92 | !!---------------------------------------------------------------------- |
---|
93 | REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: pt_su ! ice surface temperature (Kelvin) |
---|
94 | REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: ph_ice ! sea-ice thickness |
---|
95 | REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: ph_snw ! snow depth |
---|
96 | LOGICAL , INTENT(in ) :: ld_pnd_alb ! effect of melt ponds on albedo |
---|
97 | REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: pafrac_pnd ! melt pond relative fraction (per unit ice area) |
---|
98 | REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: ph_pnd ! melt pond depth |
---|
99 | REAL(wp), INTENT(in ), DIMENSION(:,:,:) :: plh_pnd ! melt pond lid thickness |
---|
100 | REAL(wp), INTENT( out), DIMENSION(:,:,:) :: palb_cs ! albedo of ice under clear sky |
---|
101 | REAL(wp), INTENT( out), DIMENSION(:,:,:) :: palb_os ! albedo of ice under overcast sky |
---|
102 | ! |
---|
103 | REAL(wp), PARAMETER :: pnd_lid_max = 0.015 ! pond lid thickness above which the ponds disappear from the albedo calculation |
---|
104 | REAL(wp), PARAMETER :: pnd_lid_min = 0.005 ! pond lid thickness below which the full pond area is used in the albedo calculation |
---|
105 | ! Note: these two variables are mirrored in sbccpl.F90 (maybe put them in one place...) |
---|
106 | ! |
---|
107 | INTEGER :: ji, jj, jl ! dummy loop indices |
---|
108 | REAL(wp) :: z1_c1, z1_c2,z1_c3, z1_c4 ! local scalar |
---|
109 | REAL(wp) :: z1_href_pnd ! inverse of the characteristic length scale (Lecomte et al. 2015) |
---|
110 | REAL(wp) :: zalb_pnd, zafrac_pnd ! ponded sea ice albedo & relative pound fraction |
---|
111 | REAL(wp) :: zalb_ice, zafrac_ice ! bare sea ice albedo & relative ice fraction |
---|
112 | REAL(wp) :: zalb_snw, zafrac_snw ! snow-covered sea ice albedo & relative snow fraction |
---|
113 | REAL(wp) :: lfrac_pnd ! The fraction of the meltpond exposed (not inder a frozen lid) |
---|
114 | !!--------------------------------------------------------------------- |
---|
115 | ! |
---|
116 | IF( ln_timing ) CALL timing_start('icealb') |
---|
117 | ! |
---|
118 | z1_href_pnd = 1. / 0.05 |
---|
119 | z1_c1 = 1. / ( LOG(1.5) - LOG(0.05) ) |
---|
120 | z1_c2 = 1. / 0.05 |
---|
121 | z1_c3 = 1. / 0.02 |
---|
122 | z1_c4 = 1. / 0.03 |
---|
123 | ! |
---|
124 | DO jl = 1, jpl |
---|
125 | DO jj = 1, jpj |
---|
126 | DO ji = 1, jpi |
---|
127 | ! !--- Specific snow, ice and pond fractions (for now, we prevent melt ponds and snow at the same time) |
---|
128 | IF( ph_snw(ji,jj,jl) == 0._wp ) THEN |
---|
129 | zafrac_snw = 0._wp |
---|
130 | IF( ld_pnd_alb ) THEN |
---|
131 | IF ( plh_pnd(ji,jj,jl) > pnd_lid_max ) THEN |
---|
132 | lfrac_pnd = 0._wp |
---|
133 | ELSE |
---|
134 | IF ( plh_pnd(ji,jj,jl) < pnd_lid_min ) THEN |
---|
135 | lfrac_pnd = 1._wp |
---|
136 | ELSE |
---|
137 | lfrac_pnd = ( plh_pnd(ji,jj,jl) - pnd_lid_min ) / (pnd_lid_max - pnd_lid_min) |
---|
138 | END IF |
---|
139 | END IF |
---|
140 | zafrac_pnd = pafrac_pnd(ji,jj,jl) * lfrac_pnd |
---|
141 | ELSE |
---|
142 | zafrac_pnd = 0._wp |
---|
143 | ENDIF |
---|
144 | zafrac_ice = 1._wp - zafrac_pnd |
---|
145 | ELSE |
---|
146 | zafrac_snw = 1._wp ! Snow fully "shades" melt ponds and ice |
---|
147 | zafrac_pnd = 0._wp |
---|
148 | zafrac_ice = 0._wp |
---|
149 | ENDIF |
---|
150 | ! |
---|
151 | ! !--- Bare ice albedo (for hi > 150cm) |
---|
152 | IF( ld_pnd_alb ) THEN |
---|
153 | zalb_ice = rn_alb_idry |
---|
154 | ELSE |
---|
155 | IF( ph_snw(ji,jj,jl) == 0._wp .AND. pt_su(ji,jj,jl) >= rt0 ) THEN ; zalb_ice = rn_alb_imlt |
---|
156 | ELSE ; zalb_ice = rn_alb_idry ; ENDIF |
---|
157 | ENDIF |
---|
158 | ! !--- Bare ice albedo (for hi < 150cm) |
---|
159 | IF( 0.05 < ph_ice(ji,jj,jl) .AND. ph_ice(ji,jj,jl) <= 1.5 ) THEN ! 5cm < hi < 150cm |
---|
160 | zalb_ice = zalb_ice + ( 0.18 - zalb_ice ) * z1_c1 * ( LOG(1.5) - LOG(ph_ice(ji,jj,jl)) ) |
---|
161 | ELSEIF( ph_ice(ji,jj,jl) <= 0.05 ) THEN ! 0cm < hi < 5cm |
---|
162 | zalb_ice = rn_alb_oce + ( 0.18 - rn_alb_oce ) * z1_c2 * ph_ice(ji,jj,jl) |
---|
163 | ENDIF |
---|
164 | ! |
---|
165 | ! !--- Snow-covered ice albedo (freezing, melting cases) |
---|
166 | IF( pt_su(ji,jj,jl) < rt0 ) THEN |
---|
167 | zalb_snw = rn_alb_sdry - ( rn_alb_sdry - zalb_ice ) * EXP( - ph_snw(ji,jj,jl) * z1_c3 ) |
---|
168 | ELSE |
---|
169 | zalb_snw = rn_alb_smlt - ( rn_alb_smlt - zalb_ice ) * EXP( - ph_snw(ji,jj,jl) * z1_c4 ) |
---|
170 | ENDIF |
---|
171 | ! !--- Ponded ice albedo |
---|
172 | IF( ld_pnd_alb ) THEN |
---|
173 | zalb_pnd = rn_alb_dpnd - ( rn_alb_dpnd - zalb_ice ) * EXP( - ph_pnd(ji,jj,jl) * z1_href_pnd ) |
---|
174 | ELSE |
---|
175 | zalb_pnd = rn_alb_dpnd |
---|
176 | ENDIF |
---|
177 | ! !--- Surface albedo is weighted mean of snow, ponds and bare ice contributions |
---|
178 | palb_os(ji,jj,jl) = ( zafrac_snw * zalb_snw + zafrac_pnd * zalb_pnd + zafrac_ice * zalb_ice ) * tmask(ji,jj,1) |
---|
179 | ! |
---|
180 | palb_cs(ji,jj,jl) = palb_os(ji,jj,jl) & |
---|
181 | & - ( - 0.1010 * palb_os(ji,jj,jl) * palb_os(ji,jj,jl) & |
---|
182 | & + 0.1933 * palb_os(ji,jj,jl) - 0.0148 ) * tmask(ji,jj,1) |
---|
183 | ! |
---|
184 | END DO |
---|
185 | END DO |
---|
186 | END DO |
---|
187 | ! |
---|
188 | ! |
---|
189 | IF( ln_timing ) CALL timing_stop('icealb') |
---|
190 | ! |
---|
191 | END SUBROUTINE ice_alb |
---|
192 | |
---|
193 | |
---|
194 | SUBROUTINE ice_alb_init |
---|
195 | !!---------------------------------------------------------------------- |
---|
196 | !! *** ROUTINE alb_init *** |
---|
197 | !! |
---|
198 | !! ** Purpose : initializations for the albedo parameters |
---|
199 | !! |
---|
200 | !! ** Method : Read the namelist namalb |
---|
201 | !!---------------------------------------------------------------------- |
---|
202 | INTEGER :: ios ! Local integer output status for namelist read |
---|
203 | !! |
---|
204 | NAMELIST/namalb/ rn_alb_sdry, rn_alb_smlt, rn_alb_idry, rn_alb_imlt, rn_alb_dpnd |
---|
205 | !!---------------------------------------------------------------------- |
---|
206 | ! |
---|
207 | REWIND( numnam_ice_ref ) ! Namelist namalb in reference namelist : Albedo parameters |
---|
208 | READ ( numnam_ice_ref, namalb, IOSTAT = ios, ERR = 901) |
---|
209 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namalb in reference namelist', lwp ) |
---|
210 | REWIND( numnam_ice_cfg ) ! Namelist namalb in configuration namelist : Albedo parameters |
---|
211 | READ ( numnam_ice_cfg, namalb, IOSTAT = ios, ERR = 902 ) |
---|
212 | 902 IF( ios > 0 ) CALL ctl_nam ( ios , 'namalb in configuration namelist', lwp ) |
---|
213 | IF(lwm) WRITE( numoni, namalb ) |
---|
214 | ! |
---|
215 | IF(lwp) THEN ! Control print |
---|
216 | WRITE(numout,*) |
---|
217 | WRITE(numout,*) 'ice_alb_init: set albedo parameters' |
---|
218 | WRITE(numout,*) '~~~~~~~~~~~~' |
---|
219 | WRITE(numout,*) ' Namelist namalb:' |
---|
220 | WRITE(numout,*) ' albedo of dry snow rn_alb_sdry = ', rn_alb_sdry |
---|
221 | WRITE(numout,*) ' albedo of melting snow rn_alb_smlt = ', rn_alb_smlt |
---|
222 | WRITE(numout,*) ' albedo of dry ice rn_alb_idry = ', rn_alb_idry |
---|
223 | WRITE(numout,*) ' albedo of bare puddled ice rn_alb_imlt = ', rn_alb_imlt |
---|
224 | WRITE(numout,*) ' albedo of ponded ice rn_alb_dpnd = ', rn_alb_dpnd |
---|
225 | ENDIF |
---|
226 | ! |
---|
227 | END SUBROUTINE ice_alb_init |
---|
228 | |
---|
229 | #else |
---|
230 | !!---------------------------------------------------------------------- |
---|
231 | !! Default option Dummy module NO SI3 sea-ice model |
---|
232 | !!---------------------------------------------------------------------- |
---|
233 | #endif |
---|
234 | |
---|
235 | !!====================================================================== |
---|
236 | END MODULE icealb |
---|