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 trunk/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: trunk/NEMOGCM/NEMO/OPA_SRC/SBC/albedo.F90 @ 7698

Last change on this file since 7698 was 7698, checked in by mocavero, 7 years ago

update trunk with OpenMP parallelization

  • Property svn:keywords set to Id
File size: 15.3 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_albice
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)            ::   ralb_im, ralb_sf, ralb_sm, ralb_if
92      REAL(wp)            ::   zswitch, z1_c1, z1_c2
93      REAL(wp)                            ::   zalb_sm, zalb_sf, zalb_st ! albedo of snow melting, freezing, total
94      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zalb, zalb_it             ! intermediate variable & albedo of ice (snow free)
95      !!---------------------------------------------------------------------
96
97      ijpl = SIZE( pt_ice, 3 )                     ! number of ice categories
98     
99      CALL wrk_alloc( jpi,jpj,ijpl, zalb, zalb_it )
100
101      IF( albd_init == 0 )   CALL albedo_init      ! initialization
102
103     
104      SELECT CASE ( nn_ice_alb )
105
106      !------------------------------------------
107      !  Shine and Henderson-Sellers (1985)
108      !------------------------------------------
109      CASE( 0 )
110       
111         ralb_sf = 0.80       ! dry snow
112         ralb_sm = 0.65       ! melting snow
113         ralb_if = 0.72       ! bare frozen ice
114         ralb_im = rn_albice  ! bare puddled ice
115         
116         !  Computation of ice albedo (free of snow)
117!$OMP PARALLEL DO schedule(static) private(jl,jj,ji)
118         DO jl = 1, ijpl
119            DO jj = 1, jpj
120               DO ji = 1, jpi
121                  IF ( ph_snw(ji,jj,jl) == 0._wp .AND. pt_ice(ji,jj,jl) >= rt0_ice ) THEN
122                     zalb(ji,jj,jl) = ralb_im
123                  ELSE
124                     zalb(ji,jj,jl) = ralb_if
125                  END IF
126               END DO
127            END DO
128         END DO
129     
130         WHERE     ( 1.5  < ph_ice                     )  ;  zalb_it = zalb
131         ELSE WHERE( 1.0  < ph_ice .AND. ph_ice <= 1.5 )  ;  zalb_it = 0.472  + 2.0 * ( zalb - 0.472 ) * ( ph_ice - 1.0 )
132         ELSE WHERE( 0.05 < ph_ice .AND. ph_ice <= 1.0 )  ;  zalb_it = 0.2467 + 0.7049 * ph_ice              &
133            &                                                                 - 0.8608 * ph_ice * ph_ice     &
134            &                                                                 + 0.3812 * ph_ice * ph_ice * ph_ice
135         ELSE WHERE                                       ;  zalb_it = 0.1    + 3.6    * ph_ice
136         END WHERE
137!$OMP PARALLEL
138!$OMP DO schedule(static) private(jl, jj, ji,zswitch,zalb_sf,zalb_sm,zalb_st)
139         DO jl = 1, ijpl
140            DO jj = 1, jpj
141               DO ji = 1, jpi
142                  ! freezing snow
143                  ! no effect of underlying ice layer IF snow thickness > c1. Albedo does not depend on snow thick if > c2
144                  !                                        !  freezing snow       
145                  zswitch   = 1._wp - MAX( 0._wp , SIGN( 1._wp , - ( ph_snw(ji,jj,jl) - c1 ) ) )
146                  zalb_sf   = ( 1._wp - zswitch ) * (  zalb_it(ji,jj,jl)  &
147                     &                           + ph_snw(ji,jj,jl) * ( ralb_sf - zalb_it(ji,jj,jl) ) / c1  )   &
148                     &        +         zswitch   * ralb_sf 
149
150                  ! melting snow
151                  ! no effect of underlying ice layer. Albedo does not depend on snow thick IF > c2
152                  zswitch   = MAX( 0._wp , SIGN( 1._wp , ph_snw(ji,jj,jl) - c2 ) )
153                  zalb_sm = ( 1._wp - zswitch ) * ( ralb_im + ph_snw(ji,jj,jl) * ( ralb_sm - ralb_im ) / c2 )   &
154                      &     +         zswitch   *   ralb_sm 
155                  !
156                  ! snow albedo
157                  zswitch  =  MAX( 0._wp , SIGN( 1._wp , pt_ice(ji,jj,jl) - rt0_snow ) )   
158                  zalb_st  =  zswitch * zalb_sm + ( 1._wp - zswitch ) * zalb_sf
159               
160                  ! Ice/snow albedo
161                  zswitch   = 1._wp - MAX( 0._wp , SIGN( 1._wp , - ph_snw(ji,jj,jl) ) )
162                  pa_ice_cs(ji,jj,jl) =  zswitch * zalb_st + ( 1._wp - zswitch ) * zalb_it(ji,jj,jl)
163                  !
164               END DO
165            END DO
166         END DO
167
168!$OMP DO schedule(static) private(jl, jj, ji)     
169         DO jl = 1, ijpl
170            DO jj = 1, jpj
171               DO ji = 1, jpi
172                  pa_ice_os(ji,jj,jl) = pa_ice_cs(ji,jj,jl) + rcloud       ! Oberhuber correction for overcast sky
173               END DO
174            END DO
175         END DO
176!$OMP END PARALLEL
177
178      !------------------------------------------
179      !  New parameterization (2016)
180      !------------------------------------------
181      CASE( 1 ) 
182
183         ralb_im = rn_albice  ! bare puddled ice
184! compilation of values from literature
185         ralb_sf = 0.85      ! dry snow
186         ralb_sm = 0.75      ! melting snow
187         ralb_if = 0.60      ! bare frozen ice
188! Perovich et al 2002 (Sheba) => the only dataset for which all types of ice/snow were retrieved
189!         ralb_sf = 0.85       ! dry snow
190!         ralb_sm = 0.72       ! melting snow
191!         ralb_if = 0.65       ! bare frozen ice
192! Brandt et al 2005 (East Antarctica)
193!         ralb_sf = 0.87      ! dry snow
194!         ralb_sm = 0.82      ! melting snow
195!         ralb_if = 0.54      ! bare frozen ice
196!
197         !  Computation of ice albedo (free of snow)
198         z1_c1 = 1. / ( LOG(1.5) - LOG(0.05) ) 
199         z1_c2 = 1. / 0.05
200         WHERE     ( ph_snw == 0._wp .AND. pt_ice >= rt0_ice )   ;   zalb = ralb_im
201         ELSE WHERE                                              ;   zalb = ralb_if
202         END  WHERE
203         
204         WHERE     ( 1.5  < ph_ice                     )  ;  zalb_it = zalb
205         ELSE WHERE( 0.05 < ph_ice .AND. ph_ice <= 1.5 )  ;  zalb_it = zalb     + ( 0.18 - zalb     ) * z1_c1 *  &
206            &                                                                     ( LOG(1.5) - LOG(ph_ice) )
207         ELSE WHERE                                       ;  zalb_it = ralb_oce + ( 0.18 - ralb_oce ) * z1_c2 * ph_ice
208         END WHERE
209
210         z1_c1 = 1. / 0.02
211         z1_c2 = 1. / 0.03
212         !  Computation of the snow/ice albedo
213!$OMP PARALLEL DO schedule(static) private(jl, jj, ji,zswitch,zalb_sf,zalb_sm,zalb_st)     
214         DO jl = 1, ijpl
215            DO jj = 1, jpj
216               DO ji = 1, jpi
217                  zalb_sf = ralb_sf - ( ralb_sf - zalb_it(ji,jj,jl)) * EXP( - ph_snw(ji,jj,jl) * z1_c1 );
218                  zalb_sm = ralb_sm - ( ralb_sm - zalb_it(ji,jj,jl)) * EXP( - ph_snw(ji,jj,jl) * z1_c2 );
219
220                   ! snow albedo
221                  zswitch = MAX( 0._wp , SIGN( 1._wp , pt_ice(ji,jj,jl) - rt0_snow ) )   
222                  zalb_st = zswitch * zalb_sm + ( 1._wp - zswitch ) * zalb_sf
223
224                  ! Ice/snow albedo   
225                  zswitch             = MAX( 0._wp , SIGN( 1._wp , - ph_snw(ji,jj,jl) ) )
226                  pa_ice_os(ji,jj,jl) = ( 1._wp - zswitch ) * zalb_st + zswitch *  zalb_it(ji,jj,jl)
227
228              END DO
229            END DO
230         END DO
231         ! Effect of the clouds (2d order polynomial)
232         pa_ice_cs = pa_ice_os - ( - 0.1010 * pa_ice_os * pa_ice_os + 0.1933 * pa_ice_os - 0.0148 ); 
233
234      END SELECT
235     
236      CALL wrk_dealloc( jpi,jpj,ijpl, zalb, zalb_it )
237      !
238   END SUBROUTINE albedo_ice
239
240
241   SUBROUTINE albedo_oce( pa_oce_os , pa_oce_cs )
242      !!----------------------------------------------------------------------
243      !!               ***  ROUTINE albedo_oce  ***
244      !!
245      !! ** Purpose :   Computation of the albedo of the ocean
246      !!----------------------------------------------------------------------
247      REAL(wp), DIMENSION(:,:), INTENT(out) ::   pa_oce_os   !  albedo of ocean under overcast sky
248      REAL(wp), DIMENSION(:,:), INTENT(out) ::   pa_oce_cs   !  albedo of ocean under clear sky
249      !!
250      REAL(wp) :: zcoef 
251      INTEGER  ::   ji, jj                                   ! dummy loop indices
252      !!----------------------------------------------------------------------
253      !
254      zcoef = 0.05 / ( 1.1 * rmue**1.4 + 0.15 )   ! Parameterization of Briegled and Ramanathan, 1982
255!$OMP PARALLEL DO schedule(static) private(jj, ji)     
256      DO jj = 1, jpj
257         DO ji = 1, jpi
258            pa_oce_cs(ji,jj) = zcoef 
259            pa_oce_os(ji,jj) = 0.06                       ! Parameterization of Kondratyev, 1969 and Payne, 1972
260         END DO
261      END DO
262      !
263   END SUBROUTINE albedo_oce
264
265
266   SUBROUTINE albedo_init
267      !!----------------------------------------------------------------------
268      !!                 ***  ROUTINE albedo_init  ***
269      !!
270      !! ** Purpose :   initializations for the albedo parameters
271      !!
272      !! ** Method  :   Read the namelist namsbc_alb
273      !!----------------------------------------------------------------------
274      INTEGER  ::   ios                 ! Local integer output status for namelist read
275      NAMELIST/namsbc_alb/ nn_ice_alb, rn_albice 
276      !!----------------------------------------------------------------------
277      !
278      albd_init = 1                     ! indicate that the initialization has been done
279      !
280      REWIND( numnam_ref )              ! Namelist namsbc_alb in reference namelist : Albedo parameters
281      READ  ( numnam_ref, namsbc_alb, IOSTAT = ios, ERR = 901)
282901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_alb in reference namelist', lwp )
283
284      REWIND( numnam_cfg )              ! Namelist namsbc_alb in configuration namelist : Albedo parameters
285      READ  ( numnam_cfg, namsbc_alb, IOSTAT = ios, ERR = 902 )
286902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namsbc_alb in configuration namelist', lwp )
287      IF(lwm) WRITE ( numond, namsbc_alb )
288      !
289      IF(lwp) THEN                      ! Control print
290         WRITE(numout,*)
291         WRITE(numout,*) 'albedo : set albedo parameters'
292         WRITE(numout,*) '~~~~~~~'
293         WRITE(numout,*) '   Namelist namsbc_alb : albedo '
294         WRITE(numout,*) '      choose the albedo parameterization                  nn_ice_alb = ', nn_ice_alb
295         WRITE(numout,*) '      albedo of bare puddled ice                          rn_albice  = ', rn_albice
296      ENDIF
297      !
298   END SUBROUTINE albedo_init
299
300   !!======================================================================
301END MODULE albedo
Note: See TracBrowser for help on using the repository browser.