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

source: trunk/NEMO/OPA_SRC/SBC/albedo.F90 @ 1208

Last change on this file since 1208 was 1208, checked in by rblod, 16 years ago

Correct syntax error with the conv, thanks to Steven, see ticket #261

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 10.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  !  01-04  (LIM 1.0)
7   !!            8.5  !  03-07  (C. Ethe, G. Madec)  Optimization (old name:shine)
8   !!            9.0  !  04-11  (C. Talandier)  add albedo_init
9   !!             -   !  01-06  (M. Vancoppenolle) LIM 3.0
10   !!             -   !  06-08  (G. Madec)  cleaning for surface module
11   !!----------------------------------------------------------------------
12   !!   albedo_ice  : albedo for   ice (clear and overcast skies)
13   !!   albedo_oce  : albedo for ocean (clear and overcast skies)
14   !!   albedo_init : initialisation of albedo computation
15   !!----------------------------------------------------------------------
16   USE phycst          ! physical constants
17   USE in_out_manager  ! I/O manager
18
19   IMPLICIT NONE
20   PRIVATE
21
22   PUBLIC albedo_ice   ! routine called sbcice_lim.F90
23   PUBLIC albedo_oce   ! routine called by ???
24
25   INTEGER  ::   albd_init = 0      !: control flag for initialization
26   REAL(wp) ::   zzero     = 0.e0   ! constant values
27   REAL(wp) ::   zone      = 1.e0   !    "       "
28
29   REAL(wp) ::   c1     = 0.05    ! constants values
30   REAL(wp) ::   c2     = 0.10    !    "        "
31   REAL(wp) ::   rmue   = 0.40    !  cosine of local solar altitude
32
33   !!* namelist namalb
34   REAL(wp) ::   &
35      cgren  = 0.06  ,     &   !  cloudiness effect on snow or ice albedo (Grenfell & Perovich, 1984)
36#if defined key_lim3
37      albice = 0.53  ,     &   !  albedo of melting ice in the arctic and antarctic (Shine & Hendersson-Sellers)
38#else
39      albice = 0.50  ,     &   !  albedo of melting ice in the arctic and antarctic (Shine & Hendersson-Sellers)
40#endif
41      alphd  = 0.80  ,     &   !  coefficients for linear interpolation used to compute
42      alphdi = 0.72  ,     &   !  albedo between two extremes values (Pyane, 1972)
43      alphc  = 0.65 
44
45   !!----------------------------------------------------------------------
46   !!   OPA 9.0 , LOCEAN-IPSL (2006)
47   !! $Id$
48   !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt)
49   !!----------------------------------------------------------------------
50
51CONTAINS
52
53   SUBROUTINE albedo_ice( pt_ice, ph_ice, ph_snw, pa_ice_cs, pa_ice_os )
54      !!----------------------------------------------------------------------
55      !!               ***  ROUTINE albedo_ice  ***
56      !!         
57      !! ** Purpose :   Computation of the albedo of the snow/ice system
58      !!                as well as the ocean one
59      !!       
60      !! ** Method  : - Computation of the albedo of snow or ice (choose the
61      !!                rignt one by a large number of tests
62      !!              - Computation of the albedo of the ocean
63      !!
64      !! References :   Shine and Hendersson-Sellers 1985, JGR, 90(D1), 2243-2250.
65      !!----------------------------------------------------------------------
66      REAL(wp), INTENT(in   ), DIMENSION(:,:,:) ::   pt_ice      !  ice surface temperature
67      REAL(wp), INTENT(in   ), DIMENSION(:,:,:) ::   ph_ice      !  sea-ice thickness
68      REAL(wp), INTENT(in   ), DIMENSION(:,:,:) ::   ph_snw      !  snow thickness
69      REAL(wp), INTENT(  out), DIMENSION(:,:,:) ::   pa_ice_cs   !  albedo of ice under clear    sky
70      REAL(wp), INTENT(  out), DIMENSION(:,:,:) ::   pa_ice_os   !  albedo of ice under overcast sky
71      !!
72      INTEGER  ::   ji, jj, jl    ! dummy loop indices
73      INTEGER  ::   ijpl          ! number of ice categories (3rd dim of ice input arrays)
74      REAL(wp) ::   zalbpsnm      ! albedo of ice under clear sky when snow is melting
75      REAL(wp) ::   zalbpsnf      ! albedo of ice under clear sky when snow is freezing
76      REAL(wp) ::   zalbpsn       ! albedo of snow/ice system when ice is coverd by snow
77      REAL(wp) ::   zalbpic       ! albedo of snow/ice system when ice is free of snow
78      REAL(wp) ::   zithsn        ! = 1 for hsn >= 0 ( ice is cov. by snow ) ; = 0 otherwise (ice is free of snow)
79      REAL(wp) ::   zitmlsn       ! = 1 freezinz snow (pt_ice >=rt0_snow) ; = 0 melting snow (pt_ice<rt0_snow)
80      REAL(wp) ::   zihsc1        ! = 1 hsn <= c1 ; = 0 hsn > c1
81      REAL(wp) ::   zihsc2        ! = 1 hsn >= c2 ; = 0 hsn < c2
82      !!
83      LOGICAL , DIMENSION(jpi,jpj,SIZE(pt_ice,3)) ::   llmask
84      REAL(wp), DIMENSION(jpi,jpj,SIZE(pt_ice,3)) ::   zalbfz    ! = alphdi for freezing ice ; = albice for melting ice
85      REAL(wp), DIMENSION(jpi,jpj,SIZE(pt_ice,3)) ::   zficeth   !  function of ice thickness
86      !!---------------------------------------------------------------------
87     
88      ijpl = SIZE( pt_ice, 3 )                     ! number of ice categories
89
90      IF( albd_init == 0 )   CALL albedo_init      ! initialization
91
92      !---------------------------
93      !  Computation of  zficeth
94      !---------------------------
95      llmask = ( ph_snw == 0.e0 ) .AND. ( pt_ice >= rt0_ice )
96      ! ice free of snow and melts
97      WHERE( llmask )   
98         zalbfz = albice
99      ELSEWHERE     
100         zalbfz = alphdi
101      END WHERE
102
103      DO jl = 1, ijpl
104         DO jj = 1, jpj
105            DO ji = 1, jpi
106               IF( ph_ice(ji,jj,jl) > 1.5 ) THEN
107                  zficeth(ji,jj,jl) = zalbfz(ji,jj,jl)
108               ELSEIF( ph_ice(ji,jj,jl) > 1.0  .AND. ph_ice(ji,jj,jl) <= 1.5 ) THEN
109                  zficeth(ji,jj,jl) = 0.472 + 2.0 * ( zalbfz(ji,jj,jl) - 0.472 ) * ( ph_ice(ji,jj,jl) - 1.0 )
110               ELSEIF( ph_ice(ji,jj,jl) > 0.05 .AND. ph_ice(ji,jj,jl) <= 1.0 ) THEN
111                  zficeth(ji,jj,jl) = 0.2467 + 0.7049 * ph_ice(ji,jj,jl)                               &
112                     &                    - 0.8608 * ph_ice(ji,jj,jl) * ph_ice(ji,jj,jl)                 &
113                     &                    + 0.3812 * ph_ice(ji,jj,jl) * ph_ice(ji,jj,jl) * ph_ice (ji,jj,jl)
114               ELSE
115                  zficeth(ji,jj,jl) = 0.1 + 3.6 * ph_ice(ji,jj,jl) 
116               ENDIF
117            END DO
118         END DO
119      END DO
120     
121      !-----------------------------------------------
122      !    Computation of the snow/ice albedo system
123      !-------------------------- ---------------------
124     
125      !    Albedo of snow-ice for clear sky.
126      !-----------------------------------------------   
127      DO jl = 1, ijpl
128         DO jj = 1, jpj
129            DO ji = 1, jpi
130               !  Case of ice covered by snow.             
131               !                                        !  freezing snow       
132               zihsc1   = 1.0 - MAX( zzero , SIGN( zone , - ( ph_snw(ji,jj,jl) - c1 ) ) )
133               zalbpsnf = ( 1.0 - zihsc1 ) * (  zficeth(ji,jj,jl)                                        &
134                  &                           + ph_snw(ji,jj,jl) * ( alphd - zficeth(ji,jj,jl) ) / c1  )   &
135                  &     +         zihsc1   * alphd 
136               !                                        !  melting snow               
137               zihsc2   = MAX( zzero , SIGN( zone , ph_snw(ji,jj,jl) - c2 ) )
138               zalbpsnm = ( 1.0 - zihsc2 ) * ( albice + ph_snw(ji,jj,jl) * ( alphc - albice ) / c2 )       &
139                  &     +         zihsc2   * alphc 
140               !
141               zitmlsn  =  MAX( zzero , SIGN( zone , pt_ice(ji,jj,jl) - rt0_snow ) )   
142               zalbpsn  =  zitmlsn * zalbpsnm + ( 1.0 - zitmlsn ) * zalbpsnf
143           
144               !  Case of ice free of snow.
145               zalbpic  = zficeth(ji,jj,jl) 
146           
147               ! albedo of the system   
148               zithsn   = 1.0 - MAX( zzero , SIGN( zone , - ph_snw(ji,jj,jl) ) )
149               pa_ice_cs(ji,jj,jl) =  zithsn * zalbpsn + ( 1.0 - zithsn ) *  zalbpic
150            END DO
151         END DO
152      END DO
153     
154      !    Albedo of snow-ice for overcast sky.
155      !---------------------------------------------- 
156      pa_ice_os(:,:,:) = pa_ice_cs(:,:,:) + cgren       ! Oberhuber correction
157      !
158   END SUBROUTINE albedo_ice
159
160
161   SUBROUTINE albedo_oce( pa_oce_os , pa_oce_cs )
162      !!----------------------------------------------------------------------
163      !!               ***  ROUTINE albedo_oce  ***
164      !!
165      !! ** Purpose :   Computation of the albedo of the ocean
166      !!
167      !! ** Method  :   ....
168      !!----------------------------------------------------------------------
169      REAL(wp), DIMENSION(jpi,jpj), INTENT(out) ::   pa_oce_os   !  albedo of ocean under overcast sky
170      REAL(wp), DIMENSION(jpi,jpj), INTENT(out) ::   pa_oce_cs   !  albedo of ocean under clear sky
171      !!
172      REAL(wp) ::   zcoef   ! temporary scalar
173      !!----------------------------------------------------------------------
174      !
175      zcoef = 0.05 / ( 1.1 * rmue**1.4 + 0.15 )      ! Parameterization of Briegled and Ramanathan, 1982
176      pa_oce_cs(:,:) = zcoef               
177      pa_oce_os(:,:)  = 0.06                         ! Parameterization of Kondratyev, 1969 and Payne, 1972
178      !
179   END SUBROUTINE albedo_oce
180
181
182   SUBROUTINE albedo_init
183      !!----------------------------------------------------------------------
184      !!                 ***  ROUTINE albedo_init  ***
185      !!
186      !! ** Purpose :   initializations for the albedo parameters
187      !!
188      !! ** Method  :   Read the namelist namalb
189      !!----------------------------------------------------------------------
190      NAMELIST/namalb/ cgren, albice, alphd, alphdi, alphc
191      !!----------------------------------------------------------------------
192
193      ! set the initialization flag to 1
194      albd_init = 1           ! indicate that the initialization has been done
195
196      ! Read Namelist namalb : albedo parameters
197      REWIND( numnam )
198      READ  ( numnam, namalb )
199
200      IF(lwp) THEN               ! Control print
201         WRITE(numout,*)
202         WRITE(numout,*) 'albedo_init : set albedo parameters from namelist namalb'
203         WRITE(numout,*) '~~~~~~~~~~~'
204         WRITE(numout,*) '             correction for snow and ice albedo                    cgren  = ', cgren
205         WRITE(numout,*) '             albedo of melting ice in the arctic and antarctic     albice = ', albice
206         WRITE(numout,*) '             coefficients for linear                               alphd  = ', alphd
207         WRITE(numout,*) '             interpolation used to compute albedo                  alphdi = ', alphdi
208         WRITE(numout,*) '             between two extremes values (Pyane, 1972)             alphc  = ', alphc
209      ENDIF
210      !
211   END SUBROUTINE albedo_init
212
213   !!======================================================================
214END MODULE albedo
Note: See TracBrowser for help on using the repository browser.