New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
albedo.F90 in branches/UKMO/dev_r8183_GC_couple_pkg/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/UKMO/dev_r8183_GC_couple_pkg/NEMOGCM/NEMO/OPA_SRC/SBC/albedo.F90 @ 8730

Last change on this file since 8730 was 8730, checked in by dancopsey, 6 years ago

Cleared out SVN keywords.

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