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 @ 2715

Last change on this file since 2715 was 2715, checked in by rblod, 13 years ago

First attempt to put dynamic allocation on the trunk

  • Property svn:keywords set to Id
File size: 11.9 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   !!----------------------------------------------------------------------
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
22   IMPLICIT NONE
23   PRIVATE
24
25   PUBLIC   albedo_ice   ! routine called sbcice_lim.F90
26   PUBLIC   albedo_oce   ! routine called by ???
27
28   INTEGER  ::   albd_init = 0      !: control flag for initialization
29   REAL(wp) ::   zzero     = 0.e0   ! constant values
30   REAL(wp) ::   zone      = 1.e0   !    "       "
31
32   REAL(wp) ::   c1     = 0.05    ! constants values
33   REAL(wp) ::   c2     = 0.10    !    "        "
34   REAL(wp) ::   rmue   = 0.40    !  cosine of local solar altitude
35
36   !                               !!* namelist namsbc_alb
37   REAL(wp) ::   rn_cloud  = 0.06   !  cloudiness effect on snow or ice albedo (Grenfell & Perovich, 1984)
38#if defined key_lim3
39   REAL(wp) ::   rn_albice = 0.53   !  albedo of melting ice in the arctic and antarctic (Shine & Hendersson-Sellers)
40#else
41   REAL(wp) ::   rn_albice = 0.50   !  albedo of melting ice in the arctic and antarctic (Shine & Hendersson-Sellers)
42#endif
43   REAL(wp) ::   rn_alphd  = 0.80   !  coefficients for linear interpolation used to compute
44   REAL(wp) ::   rn_alphdi = 0.72   !  albedo between two extremes values (Pyane, 1972)
45   REAL(wp) ::   rn_alphc  = 0.65   !
46
47   !!----------------------------------------------------------------------
48   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
49   !! $Id$
50   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
51   !!----------------------------------------------------------------------
52CONTAINS
53
54   SUBROUTINE albedo_ice( pt_ice, ph_ice, ph_snw, pa_ice_cs, pa_ice_os )
55      !!----------------------------------------------------------------------
56      !!               ***  ROUTINE albedo_ice  ***
57      !!         
58      !! ** Purpose :   Computation of the albedo of the snow/ice system
59      !!                as well as the ocean one
60      !!       
61      !! ** Method  : - Computation of the albedo of snow or ice (choose the
62      !!                rignt one by a large number of tests
63      !!              - Computation of the albedo of the ocean
64      !!
65      !! References :   Shine and Hendersson-Sellers 1985, JGR, 90(D1), 2243-2250.
66      !!----------------------------------------------------------------------
67      USE wrk_nemo, ONLY:   wrk_in_use, wrk_not_released
68      USE wrk_nemo, ONLY:   wrk_3d_6 , wrk_3d_7    ! 3D workspace
69      !!
70      REAL(wp), INTENT(in   ), DIMENSION(:,:,:) ::   pt_ice      !  ice surface temperature (Kelvin)
71      REAL(wp), INTENT(in   ), DIMENSION(:,:,:) ::   ph_ice      !  sea-ice thickness
72      REAL(wp), INTENT(in   ), DIMENSION(:,:,:) ::   ph_snw      !  snow thickness
73      REAL(wp), INTENT(  out), DIMENSION(:,:,:) ::   pa_ice_cs   !  albedo of ice under clear    sky
74      REAL(wp), INTENT(  out), DIMENSION(:,:,:) ::   pa_ice_os   !  albedo of ice under overcast sky
75      !!
76      INTEGER  ::   ji, jj, jl    ! dummy loop indices
77      INTEGER  ::   ijpl          ! number of ice categories (3rd dim of ice input arrays)
78      REAL(wp) ::   zalbpsnm      ! albedo of ice under clear sky when snow is melting
79      REAL(wp) ::   zalbpsnf      ! albedo of ice under clear sky when snow is freezing
80      REAL(wp) ::   zalbpsn       ! albedo of snow/ice system when ice is coverd by snow
81      REAL(wp) ::   zalbpic       ! albedo of snow/ice system when ice is free of snow
82      REAL(wp) ::   zithsn        ! = 1 for hsn >= 0 ( ice is cov. by snow ) ; = 0 otherwise (ice is free of snow)
83      REAL(wp) ::   zitmlsn       ! = 1 freezinz snow (pt_ice >=rt0_snow) ; = 0 melting snow (pt_ice<rt0_snow)
84      REAL(wp) ::   zihsc1        ! = 1 hsn <= c1 ; = 0 hsn > c1
85      REAL(wp) ::   zihsc2        ! = 1 hsn >= c2 ; = 0 hsn < c2
86      !!
87      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zalbfz    ! = rn_alphdi for freezing ice ; = rn_albice for melting ice
88      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zficeth   !  function of ice thickness
89      !!---------------------------------------------------------------------
90     
91      ijpl = SIZE( pt_ice, 3 )                     ! number of ice categories
92
93      IF( wrk_in_use(3, 6,7) ) THEN
94         CALL ctl_stop('albedo_ice: requested workspace arrays are unavailable')   ;   RETURN
95      ENDIF
96      ! Associate pointers with sub-arrays of workspace arrays
97      zalbfz  =>   wrk_3d_6(:,:,1:ijpl)
98      zficeth =>   wrk_3d_7(:,:,1:ijpl)
99
100      IF( albd_init == 0 )   CALL albedo_init      ! initialization
101
102      !---------------------------
103      !  Computation of  zficeth
104      !---------------------------
105      ! ice free of snow and melts
106      WHERE     ( ph_snw == 0._wp .AND. pt_ice >= rt0_ice )   ;   zalbfz(:,:,:) = rn_albice
107      ELSE WHERE                                              ;   zalbfz(:,:,:) = rn_alphdi
108      END  WHERE
109
110      WHERE     ( 1.5  < ph_ice                     )  ;  zficeth = zalbfz
111      ELSE WHERE( 1.0  < ph_ice .AND. ph_ice <= 1.5 )  ;  zficeth = 0.472  + 2.0 * ( zalbfz - 0.472 ) * ( ph_ice - 1.0 )
112      ELSE WHERE( 0.05 < ph_ice .AND. ph_ice <= 1.0 )  ;  zficeth = 0.2467 + 0.7049 * ph_ice              &
113         &                                                                 - 0.8608 * ph_ice * ph_ice     &
114         &                                                                 + 0.3812 * ph_ice * ph_ice * ph_ice
115      ELSE WHERE                                       ;  zficeth = 0.1    + 3.6    * ph_ice
116      END WHERE
117
118!!gm old code
119!      DO jl = 1, ijpl
120!         DO jj = 1, jpj
121!            DO ji = 1, jpi
122!               IF( ph_ice(ji,jj,jl) > 1.5 ) THEN
123!                  zficeth(ji,jj,jl) = zalbfz(ji,jj,jl)
124!               ELSEIF( ph_ice(ji,jj,jl) > 1.0  .AND. ph_ice(ji,jj,jl) <= 1.5 ) THEN
125!                  zficeth(ji,jj,jl) = 0.472 + 2.0 * ( zalbfz(ji,jj,jl) - 0.472 ) * ( ph_ice(ji,jj,jl) - 1.0 )
126!               ELSEIF( ph_ice(ji,jj,jl) > 0.05 .AND. ph_ice(ji,jj,jl) <= 1.0 ) THEN
127!                  zficeth(ji,jj,jl) = 0.2467 + 0.7049 * ph_ice(ji,jj,jl)                               &
128!                     &                    - 0.8608 * ph_ice(ji,jj,jl) * ph_ice(ji,jj,jl)                 &
129!                     &                    + 0.3812 * ph_ice(ji,jj,jl) * ph_ice(ji,jj,jl) * ph_ice (ji,jj,jl)
130!               ELSE
131!                  zficeth(ji,jj,jl) = 0.1 + 3.6 * ph_ice(ji,jj,jl)
132!               ENDIF
133!            END DO
134!         END DO
135!      END DO
136!!gm end old code
137     
138      !-----------------------------------------------
139      !    Computation of the snow/ice albedo system
140      !-------------------------- ---------------------
141     
142      !    Albedo of snow-ice for clear sky.
143      !-----------------------------------------------   
144      DO jl = 1, ijpl
145         DO jj = 1, jpj
146            DO ji = 1, jpi
147               !  Case of ice covered by snow.             
148               !                                        !  freezing snow       
149               zihsc1   = 1.0 - MAX( zzero , SIGN( zone , - ( ph_snw(ji,jj,jl) - c1 ) ) )
150               zalbpsnf = ( 1.0 - zihsc1 ) * (  zficeth(ji,jj,jl)                                             &
151                  &                           + ph_snw(ji,jj,jl) * ( rn_alphd - zficeth(ji,jj,jl) ) / c1  )   &
152                  &     +         zihsc1   * rn_alphd 
153               !                                        !  melting snow               
154               zihsc2   = MAX( zzero , SIGN( zone , ph_snw(ji,jj,jl) - c2 ) )
155               zalbpsnm = ( 1.0 - zihsc2 ) * ( rn_albice + ph_snw(ji,jj,jl) * ( rn_alphc - rn_albice ) / c2 )   &
156                  &     +         zihsc2   *   rn_alphc 
157               !
158               zitmlsn  =  MAX( zzero , SIGN( zone , pt_ice(ji,jj,jl) - rt0_snow ) )   
159               zalbpsn  =  zitmlsn * zalbpsnm + ( 1.0 - zitmlsn ) * zalbpsnf
160           
161               !  Case of ice free of snow.
162               zalbpic  = zficeth(ji,jj,jl) 
163           
164               ! albedo of the system   
165               zithsn   = 1.0 - MAX( zzero , SIGN( zone , - ph_snw(ji,jj,jl) ) )
166               pa_ice_cs(ji,jj,jl) =  zithsn * zalbpsn + ( 1.0 - zithsn ) *  zalbpic
167            END DO
168         END DO
169      END DO
170     
171      !    Albedo of snow-ice for overcast sky.
172      !---------------------------------------------- 
173      pa_ice_os(:,:,:) = pa_ice_cs(:,:,:) + rn_cloud       ! Oberhuber correction
174      !
175      IF( wrk_not_released(3, 6,7) )   CALL ctl_stop('albedo_ice: failed to release workspace arrays')
176      !
177   END SUBROUTINE albedo_ice
178
179
180   SUBROUTINE albedo_oce( pa_oce_os , pa_oce_cs )
181      !!----------------------------------------------------------------------
182      !!               ***  ROUTINE albedo_oce  ***
183      !!
184      !! ** Purpose :   Computation of the albedo of the ocean
185      !!----------------------------------------------------------------------
186      REAL(wp), DIMENSION(:,:), INTENT(out) ::   pa_oce_os   !  albedo of ocean under overcast sky
187      REAL(wp), DIMENSION(:,:), INTENT(out) ::   pa_oce_cs   !  albedo of ocean under clear sky
188      !!
189      REAL(wp) ::   zcoef   ! local scalar
190      !!----------------------------------------------------------------------
191      !
192      zcoef = 0.05 / ( 1.1 * rmue**1.4 + 0.15 )      ! Parameterization of Briegled and Ramanathan, 1982
193      pa_oce_cs(:,:) = zcoef               
194      pa_oce_os(:,:)  = 0.06                         ! Parameterization of Kondratyev, 1969 and Payne, 1972
195      !
196   END SUBROUTINE albedo_oce
197
198
199   SUBROUTINE albedo_init
200      !!----------------------------------------------------------------------
201      !!                 ***  ROUTINE albedo_init  ***
202      !!
203      !! ** Purpose :   initializations for the albedo parameters
204      !!
205      !! ** Method  :   Read the namelist namsbc_alb
206      !!----------------------------------------------------------------------
207      NAMELIST/namsbc_alb/ rn_cloud, rn_albice, rn_alphd, rn_alphdi, rn_alphc
208      !!----------------------------------------------------------------------
209      !
210      albd_init = 1                     ! indicate that the initialization has been done
211      !
212      REWIND( numnam )                  ! Read Namelist namsbc_alb : albedo parameters
213      READ  ( numnam, namsbc_alb )
214      !
215      IF(lwp) THEN                      ! Control print
216         WRITE(numout,*)
217         WRITE(numout,*) 'albedo : set albedo parameters'
218         WRITE(numout,*) '~~~~~~~'
219         WRITE(numout,*) '   Namelist namsbc_alb : albedo '
220         WRITE(numout,*) '      correction for snow and ice albedo                  rn_cloud  = ', rn_cloud
221         WRITE(numout,*) '      albedo of melting ice in the arctic and antarctic   rn_albice = ', rn_albice
222         WRITE(numout,*) '      coefficients for linear                             rn_alphd  = ', rn_alphd
223         WRITE(numout,*) '      interpolation used to compute albedo                rn_alphdi = ', rn_alphdi
224         WRITE(numout,*) '      between two extremes values (Pyane, 1972)           rn_alphc  = ', rn_alphc
225      ENDIF
226      !
227   END SUBROUTINE albedo_init
228
229   !!======================================================================
230END MODULE albedo
Note: See TracBrowser for help on using the repository browser.