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.
Changeset 6347 for branches/2016/dev_r6325_SIMPLIF_1/NEMOGCM/NEMO/OPA_SRC/SBC/albedo.F90 – NEMO

Ignore:
Timestamp:
2016-02-24T08:56:48+01:00 (8 years ago)
Author:
gm
Message:

#1683: SIMPLIF-1 : Phase with the v3.6_Stable (DOC+ZDF+traqsr+lbedo)

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2016/dev_r6325_SIMPLIF_1/NEMOGCM/NEMO/OPA_SRC/SBC/albedo.F90

    r4624 r6347  
    88   !!             -   ! 2004-11  (C. Talandier)  add albedo_init 
    99   !!             -   ! 2001-06  (M. Vancoppenolle) LIM 3.0 
    10    !!             -   ! 2006-08  (G. Madec)  cleaning for surface module 
     10   !!            3.2  ! 2006-08  (G. Madec)  cleaning for surface module 
     11   !!            3.6  ! 2016-01  (C. Rousset) new parameterization for sea ice albedo 
    1112   !!---------------------------------------------------------------------- 
    1213 
     
    1718   !!---------------------------------------------------------------------- 
    1819   USE phycst         ! physical constants 
     20   ! 
    1921   USE in_out_manager ! I/O manager 
    2022   USE lib_mpp        ! MPP library 
     
    2527   PRIVATE 
    2628 
    27    PUBLIC   albedo_ice   ! routine called sbcice_lim.F90 
    28    PUBLIC   albedo_oce   ! routine called by ??? 
    29  
    30    INTEGER  ::   albd_init = 0      !: control flag for initialization 
    31    REAL(wp) ::   zzero     = 0.e0   ! constant values 
    32    REAL(wp) ::   zone      = 1.e0   !    "       " 
    33  
    34    REAL(wp) ::   c1     = 0.05    ! constants values 
    35    REAL(wp) ::   c2     = 0.10    !    "        " 
    36    REAL(wp) ::   rmue   = 0.40    !  cosine of local solar altitude 
    37  
     29   PUBLIC   albedo_ice   ! called sbcice_lim.F90 
     30   PUBLIC   albedo_oce   ! called by sbccpl.F90 
     31 
     32   INTEGER  ::   albd_init = 0          ! control flag for initialization 
     33   
     34   REAL(wp) ::   rmue      = 0.400_wp   ! cosine of local solar altitude 
     35   REAL(wp) ::   ralb_oce  = 0.066_wp   ! ocean or lead albedo (Pegau and Paulson, Ann. Glac. 2001) 
     36   REAL(wp) ::   c1        = 0.050_wp   ! snow thickness (only for nn_ice_alb=0) 
     37   REAL(wp) ::   c2        = 0.100_wp   !  "        " 
     38   REAL(wp) ::   rcloud    = 0.060_wp   ! cloud effect on albedo (only-for nn_ice_alb=0) 
     39  
    3840   !                             !!* namelist namsbc_alb 
    39    REAL(wp) ::   rn_cloud         !  cloudiness effect on snow or ice albedo (Grenfell & Perovich, 1984) 
    40 #if defined key_lim3 
    41    REAL(wp) ::   rn_albice        !  albedo of melting ice in the arctic and antarctic (Shine & Hendersson-Sellers) 
    42 #else 
    43    REAL(wp) ::   rn_albice        !  albedo of melting ice in the arctic and antarctic (Shine & Hendersson-Sellers) 
    44 #endif 
    45    REAL(wp) ::   rn_alphd         !  coefficients for linear interpolation used to compute 
    46    REAL(wp) ::   rn_alphdi        !  albedo between two extremes values (Pyane, 1972) 
    47    REAL(wp) ::   rn_alphc         !  
     41   INTEGER  ::   nn_ice_alb       ! type of albedo parameterization  
     42   REAL(wp) ::   rn_albice        ! albedo of bare puddled ice 
    4843 
    4944   !!---------------------------------------------------------------------- 
     
    5954      !!           
    6055      !! ** Purpose :   Computation of the albedo of the snow/ice system  
    61       !!                as well as the ocean one 
    6256      !!        
    63       !! ** Method  : - Computation of the albedo of snow or ice (choose the  
    64       !!                rignt one by a large number of tests 
    65       !!              - Computation of the albedo of the ocean 
    66       !! 
    67       !! References :   Shine and Hendersson-Sellers 1985, JGR, 90(D1), 2243-2250. 
     57      !! ** Method  :   Two schemes are available (from namelist parameter nn_ice_alb) 
     58      !!                  0: the scheme is that of Shine & Henderson-Sellers (JGR 1985) for clear-skies 
     59      !!                  1: the scheme is "home made" (for cloudy skies) and based on Brandt et al. (J. Climate 2005) 
     60      !!                                                                           and Grenfell & Perovich (JGR 2004) 
     61      !!                Description of scheme 1: 
     62      !!                  1) Albedo dependency on ice thickness follows the findings from Brandt et al (2005) 
     63      !!                     which are an update of Allison et al. (JGR 1993) ; Brandt et al. 1999 
     64      !!                     0-5cm  : linear function of ice thickness 
     65      !!                     5-150cm: log    function of ice thickness 
     66      !!                     > 150cm: constant 
     67      !!                  2) Albedo dependency on snow thickness follows the findings from Grenfell & Perovich (2004) 
     68      !!                     i.e. it increases as -EXP(-snw_thick/0.02) during freezing and -EXP(-snw_thick/0.03) during melting 
     69      !!                  3) Albedo dependency on clouds is speculated from measurements of Grenfell and Perovich (2004) 
     70      !!                     i.e. cloudy-clear albedo depend on cloudy albedo following a 2d order polynomial law 
     71      !!                  4) The needed 4 parameters are: dry and melting snow, freezing ice and bare puddled ice 
     72      !! 
     73      !! ** Note    :   The parameterization from Shine & Henderson-Sellers presents several misconstructions: 
     74      !!                  1) ice albedo when ice thick. tends to 0 is different than ocean albedo 
     75      !!                  2) for small ice thick. covered with some snow (<3cm?), albedo is larger  
     76      !!                     under melting conditions than under freezing conditions 
     77      !!                  3) the evolution of ice albedo as a function of ice thickness shows   
     78      !!                     3 sharp inflexion points (at 5cm, 100cm and 150cm) that look highly unrealistic 
     79      !! 
     80      !! References :   Shine & Henderson-Sellers 1985, JGR, 90(D1), 2243-2250. 
     81      !!                Brandt et al. 2005, J. Climate, vol 18 
     82      !!                Grenfell & Perovich 2004, JGR, vol 109  
    6883      !!---------------------------------------------------------------------- 
    6984      REAL(wp), INTENT(in   ), DIMENSION(:,:,:) ::   pt_ice      !  ice surface temperature (Kelvin) 
     
    7388      REAL(wp), INTENT(  out), DIMENSION(:,:,:) ::   pa_ice_os   !  albedo of ice under overcast sky 
    7489      !! 
    75       INTEGER  ::   ji, jj, jl    ! dummy loop indices 
    76       INTEGER  ::   ijpl          ! number of ice categories (3rd dim of ice input arrays) 
    77       REAL(wp) ::   zalbpsnm      ! albedo of ice under clear sky when snow is melting 
    78       REAL(wp) ::   zalbpsnf      ! albedo of ice under clear sky when snow is freezing 
    79       REAL(wp) ::   zalbpsn       ! albedo of snow/ice system when ice is coverd by snow 
    80       REAL(wp) ::   zalbpic       ! albedo of snow/ice system when ice is free of snow 
    81       REAL(wp) ::   zithsn        ! = 1 for hsn >= 0 ( ice is cov. by snow ) ; = 0 otherwise (ice is free of snow) 
    82       REAL(wp) ::   zitmlsn       ! = 1 freezinz snow (pt_ice >=rt0_snow) ; = 0 melting snow (pt_ice<rt0_snow) 
    83       REAL(wp) ::   zihsc1        ! = 1 hsn <= c1 ; = 0 hsn > c1 
    84       REAL(wp) ::   zihsc2        ! = 1 hsn >= c2 ; = 0 hsn < c2 
    85       !! 
    86       REAL(wp), POINTER, DIMENSION(:,:,:) ::   zalbfz    ! = rn_alphdi for freezing ice ; = rn_albice for melting ice 
    87       REAL(wp), POINTER, DIMENSION(:,:,:) ::   zficeth   !  function of ice thickness 
     90      INTEGER  ::   ji, jj, jl   ! dummy loop indices 
     91      INTEGER  ::   ijpl         ! number of ice categories (3rd dim of ice input arrays) 
     92      REAL(wp) ::   ralb_im, ralb_sf, ralb_sm, ralb_if 
     93      REAL(wp) ::   zswitch, z1_c1, z1_c2 
     94      REAL(wp) ::   zalb_sm, zalb_sf, zalb_st   ! albedo of snow melting, freezing, total 
     95      REAL(wp), POINTER, DIMENSION(:,:,:) ::   zalb, zalb_it   ! intermediate variable & albedo of ice (snow free) 
    8896      !!--------------------------------------------------------------------- 
     97 
     98      ijpl = SIZE( pt_ice, 3 )                     ! number of ice categories 
    8999       
    90       ijpl = SIZE( pt_ice, 3 )                     ! number of ice categories 
    91  
    92       CALL wrk_alloc( jpi,jpj,ijpl, zalbfz, zficeth ) 
     100      CALL wrk_alloc( jpi,jpj,ijpl,   zalb, zalb_it ) 
    93101 
    94102      IF( albd_init == 0 )   CALL albedo_init      ! initialization  
    95103 
    96       !--------------------------- 
    97       !  Computation of  zficeth 
    98       !--------------------------- 
    99       ! ice free of snow and melts 
    100       WHERE     ( ph_snw == 0._wp .AND. pt_ice >= rt0_ice )   ;   zalbfz(:,:,:) = rn_albice 
    101       ELSE WHERE                                              ;   zalbfz(:,:,:) = rn_alphdi 
    102       END  WHERE 
    103  
    104       WHERE     ( 1.5  < ph_ice                     )  ;  zficeth = zalbfz 
    105       ELSE WHERE( 1.0  < ph_ice .AND. ph_ice <= 1.5 )  ;  zficeth = 0.472  + 2.0 * ( zalbfz - 0.472 ) * ( ph_ice - 1.0 ) 
    106       ELSE WHERE( 0.05 < ph_ice .AND. ph_ice <= 1.0 )  ;  zficeth = 0.2467 + 0.7049 * ph_ice              & 
    107          &                                                                 - 0.8608 * ph_ice * ph_ice     & 
    108          &                                                                 + 0.3812 * ph_ice * ph_ice * ph_ice 
    109       ELSE WHERE                                       ;  zficeth = 0.1    + 3.6    * ph_ice 
    110       END WHERE 
    111  
    112 !!gm old code 
    113 !      DO jl = 1, ijpl 
    114 !         DO jj = 1, jpj 
    115 !            DO ji = 1, jpi 
    116 !               IF( ph_ice(ji,jj,jl) > 1.5 ) THEN 
    117 !                  zficeth(ji,jj,jl) = zalbfz(ji,jj,jl) 
    118 !               ELSEIF( ph_ice(ji,jj,jl) > 1.0  .AND. ph_ice(ji,jj,jl) <= 1.5 ) THEN 
    119 !                  zficeth(ji,jj,jl) = 0.472 + 2.0 * ( zalbfz(ji,jj,jl) - 0.472 ) * ( ph_ice(ji,jj,jl) - 1.0 ) 
    120 !               ELSEIF( ph_ice(ji,jj,jl) > 0.05 .AND. ph_ice(ji,jj,jl) <= 1.0 ) THEN 
    121 !                  zficeth(ji,jj,jl) = 0.2467 + 0.7049 * ph_ice(ji,jj,jl)                               & 
    122 !                     &                    - 0.8608 * ph_ice(ji,jj,jl) * ph_ice(ji,jj,jl)                 & 
    123 !                     &                    + 0.3812 * ph_ice(ji,jj,jl) * ph_ice(ji,jj,jl) * ph_ice (ji,jj,jl) 
    124 !               ELSE 
    125 !                  zficeth(ji,jj,jl) = 0.1 + 3.6 * ph_ice(ji,jj,jl)  
    126 !               ENDIF 
    127 !            END DO 
    128 !         END DO 
    129 !      END DO 
    130 !!gm end old code 
    131104       
    132       !-----------------------------------------------  
    133       !    Computation of the snow/ice albedo system  
    134       !-------------------------- --------------------- 
     105      SELECT CASE ( nn_ice_alb )    ! Select the albedo parameterization 
     106      ! 
     107      !              !------------------------------------------ 
     108      CASE( 0 )      !  Shine and Henderson-Sellers (1985) 
     109         !           !------------------------------------------ 
     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         WHERE     ( ph_snw == 0._wp .AND. pt_ice >= rt0_ice )   ;   zalb(:,:,:) = ralb_im 
     118         ELSE WHERE                                              ;   zalb(:,:,:) = ralb_if 
     119         END  WHERE 
    135120       
    136       !    Albedo of snow-ice for clear sky. 
    137       !-----------------------------------------------     
    138       DO jl = 1, ijpl 
    139          DO jj = 1, jpj 
    140             DO ji = 1, jpi 
    141                !  Case of ice covered by snow.              
    142                !                                        !  freezing snow         
    143                zihsc1   = 1.0 - MAX( zzero , SIGN( zone , - ( ph_snw(ji,jj,jl) - c1 ) ) ) 
    144                zalbpsnf = ( 1.0 - zihsc1 ) * (  zficeth(ji,jj,jl)                                             & 
    145                   &                           + ph_snw(ji,jj,jl) * ( rn_alphd - zficeth(ji,jj,jl) ) / c1  )   & 
    146                   &     +         zihsc1   * rn_alphd   
    147                !                                        !  melting snow                 
    148                zihsc2   = MAX( zzero , SIGN( zone , ph_snw(ji,jj,jl) - c2 ) ) 
    149                zalbpsnm = ( 1.0 - zihsc2 ) * ( rn_albice + ph_snw(ji,jj,jl) * ( rn_alphc - rn_albice ) / c2 )   & 
    150                   &     +         zihsc2   *   rn_alphc  
    151                ! 
    152                zitmlsn  =  MAX( zzero , SIGN( zone , pt_ice(ji,jj,jl) - rt0_snow ) )    
    153                zalbpsn  =  zitmlsn * zalbpsnm + ( 1.0 - zitmlsn ) * zalbpsnf 
    154              
    155                !  Case of ice free of snow. 
    156                zalbpic  = zficeth(ji,jj,jl)  
    157              
    158                ! albedo of the system    
    159                zithsn   = 1.0 - MAX( zzero , SIGN( zone , - ph_snw(ji,jj,jl) ) ) 
    160                pa_ice_cs(ji,jj,jl) =  zithsn * zalbpsn + ( 1.0 - zithsn ) *  zalbpic 
     121         WHERE     ( 1.5  < ph_ice                     )  ;  zalb_it = zalb 
     122         ELSE WHERE( 1.0  < ph_ice .AND. ph_ice <= 1.5 )  ;  zalb_it = 0.472  + 2.0 * ( zalb - 0.472 ) * ( ph_ice - 1.0 ) 
     123         ELSE WHERE( 0.05 < ph_ice .AND. ph_ice <= 1.0 )  ;  zalb_it = 0.2467 + 0.7049 * ph_ice              & 
     124            &                                                                 - 0.8608 * ph_ice * ph_ice     & 
     125            &                                                                 + 0.3812 * ph_ice * ph_ice * ph_ice 
     126         ELSE WHERE                                       ;  zalb_it = 0.1    + 3.6    * ph_ice 
     127         END WHERE 
     128         ! 
     129         DO jl = 1, ijpl 
     130            DO jj = 1, jpj 
     131               DO ji = 1, jpi 
     132                  ! freezing snow 
     133                  ! no effect of underlying ice layer IF snow thickness > c1. Albedo does not depend on snow thick if > c2 
     134                  !                                        !  freezing snow         
     135                  zswitch   = 1._wp - MAX( 0._wp , SIGN( 1._wp , - ( ph_snw(ji,jj,jl) - c1 ) ) ) 
     136                  zalb_sf   = ( 1._wp - zswitch ) * (  zalb_it(ji,jj,jl)  & 
     137                     &                           + ph_snw(ji,jj,jl) * ( ralb_sf - zalb_it(ji,jj,jl) ) / c1  )   & 
     138                     &        +         zswitch   * ralb_sf   
     139                  ! 
     140                  ! melting snow 
     141                  ! no effect of underlying ice layer. Albedo does not depend on snow thick IF > c2 
     142                  zswitch   = MAX( 0._wp , SIGN( 1._wp , ph_snw(ji,jj,jl) - c2 ) ) 
     143                  zalb_sm = ( 1._wp - zswitch ) * ( ralb_im + ph_snw(ji,jj,jl) * ( ralb_sm - ralb_im ) / c2 )   & 
     144                      &     +         zswitch   *   ralb_sm  
     145                  ! 
     146                  ! snow albedo 
     147                  zswitch  =  MAX( 0._wp , SIGN( 1._wp , pt_ice(ji,jj,jl) - rt0_snow ) )    
     148                  zalb_st  =  zswitch * zalb_sm + ( 1._wp - zswitch ) * zalb_sf 
     149                  ! 
     150                  ! Ice/snow albedo 
     151                  zswitch   = 1._wp - MAX( 0._wp , SIGN( 1._wp , - ph_snw(ji,jj,jl) ) ) 
     152                  pa_ice_cs(ji,jj,jl) =  zswitch * zalb_st + ( 1._wp - zswitch ) * zalb_it(ji,jj,jl) 
     153                  ! 
     154               END DO 
    161155            END DO 
    162156         END DO 
    163       END DO 
    164        
    165       !    Albedo of snow-ice for overcast sky. 
    166       !----------------------------------------------   
    167       pa_ice_os(:,:,:) = pa_ice_cs(:,:,:) + rn_cloud       ! Oberhuber correction 
    168       ! 
    169       CALL wrk_dealloc( jpi,jpj,ijpl, zalbfz, zficeth ) 
     157         ! 
     158         pa_ice_os(:,:,:) = pa_ice_cs(:,:,:) + rcloud       ! Oberhuber correction for overcast sky 
     159         ! 
     160         !              !------------------------------------------ 
     161      CASE( 1 )         !  New parameterization (2016) 
     162         !              !------------------------------------------ 
     163         ! 
     164         ralb_im = rn_albice  ! bare puddled ice 
     165! compilation of values from literature 
     166         ralb_sf = 0.85      ! dry snow 
     167         ralb_sm = 0.75      ! melting snow 
     168         ralb_if = 0.60      ! bare frozen ice 
     169! Perovich et al 2002 (Sheba) => the only dataset for which all types of ice/snow were retrieved 
     170!         ralb_sf = 0.85       ! dry snow 
     171!         ralb_sm = 0.72       ! melting snow 
     172!         ralb_if = 0.65       ! bare frozen ice 
     173! Brandt et al 2005 (East Antarctica) 
     174!         ralb_sf = 0.87      ! dry snow 
     175!         ralb_sm = 0.82      ! melting snow 
     176!         ralb_if = 0.54      ! bare frozen ice 
     177!  
     178         !  Computation of ice albedo (free of snow) 
     179         z1_c1 = 1._wp / ( LOG(1.5_wp) - LOG(0.05_wp) )  
     180         z1_c2 = 1._wp / 0.05_wp 
     181         WHERE     ( ph_snw == 0._wp .AND. pt_ice >= rt0_ice )   ;   zalb = ralb_im 
     182         ELSE WHERE                                              ;   zalb = ralb_if 
     183         END  WHERE 
     184         ! 
     185         WHERE     ( 1.5  < ph_ice                     )  ;  zalb_it = zalb 
     186         ELSE WHERE( 0.05 < ph_ice .AND. ph_ice <= 1.5 )  ;  zalb_it = zalb     + ( 0.18_wp - zalb     ) * z1_c1 *  & 
     187            &                                                                     ( LOG(1.5_wp) - LOG(ph_ice) ) 
     188         ELSE WHERE                                       ;  zalb_it = ralb_oce + ( 0.18_wp - ralb_oce ) * z1_c2 * ph_ice 
     189         END WHERE 
     190         ! 
     191         z1_c1 = 1._wp / 0.02_wp 
     192         z1_c2 = 1._wp / 0.03_wp 
     193         !  Computation of the snow/ice albedo 
     194         DO jl = 1, ijpl 
     195            DO jj = 1, jpj 
     196               DO ji = 1, jpi 
     197                  zalb_sf = ralb_sf - ( ralb_sf - zalb_it(ji,jj,jl)) * EXP( - ph_snw(ji,jj,jl) * z1_c1 ); 
     198                  zalb_sm = ralb_sm - ( ralb_sm - zalb_it(ji,jj,jl)) * EXP( - ph_snw(ji,jj,jl) * z1_c2 ); 
     199                  ! 
     200                  ! snow albedo 
     201                  zswitch = MAX( 0._wp , SIGN( 1._wp , pt_ice(ji,jj,jl) - rt0_snow ) )    
     202                  zalb_st = zswitch * zalb_sm + ( 1._wp - zswitch ) * zalb_sf 
     203                  ! 
     204                  ! Ice/snow albedo    
     205                  zswitch             = MAX( 0._wp , SIGN( 1._wp , - ph_snw(ji,jj,jl) ) ) 
     206                  pa_ice_os(ji,jj,jl) = zswitch *  zalb_it(ji,jj,jl) + ( 1._wp - zswitch ) * zalb_st 
     207 
     208              END DO 
     209            END DO 
     210         END DO 
     211         ! Effect of the clouds (2d order polynomial) 
     212         pa_ice_cs = pa_ice_os - ( - 0.1010_wp * pa_ice_os * pa_ice_os + 0.1933_wp * pa_ice_os - 0.0148_wp )  
     213!!gm better coding ? 
     214!         pa_ice_cs = 0.1010_wp * pa_ice_os + 0.8067_wp ) * pa_ice_os + 0.0148_wp          
     215! and even better: merge this implicit do loop with the previous one 
     216!!gm end 
     217      END SELECT 
     218      ! 
     219      CALL wrk_dealloc( jpi,jpj,ijpl,   zalb, zalb_it ) 
    170220      ! 
    171221   END SUBROUTINE albedo_ice 
     
    181231      REAL(wp), DIMENSION(:,:), INTENT(out) ::   pa_oce_cs   !  albedo of ocean under clear sky 
    182232      !! 
    183       REAL(wp) ::   zcoef   ! local scalar 
    184       !!---------------------------------------------------------------------- 
    185       ! 
    186       zcoef = 0.05 / ( 1.1 * rmue**1.4 + 0.15 )      ! Parameterization of Briegled and Ramanathan, 1982  
    187       pa_oce_cs(:,:) = zcoef                
    188       pa_oce_os(:,:)  = 0.06                         ! Parameterization of Kondratyev, 1969 and Payne, 1972 
     233      REAL(wp) :: zcoef  
     234      !!---------------------------------------------------------------------- 
     235      ! 
     236      zcoef = 0.05 / ( 1.1 * rmue**1.4 + 0.15 )   ! Parameterization of Briegled and Ramanathan, 1982 
     237      pa_oce_cs(:,:) = zcoef  
     238      pa_oce_os(:,:) = 0.06                       ! Parameterization of Kondratyev, 1969 and Payne, 1972 
    189239      ! 
    190240   END SUBROUTINE albedo_oce 
     
    200250      !!---------------------------------------------------------------------- 
    201251      INTEGER  ::   ios                 ! Local integer output status for namelist read 
    202       NAMELIST/namsbc_alb/ rn_cloud, rn_albice, rn_alphd, rn_alphdi, rn_alphc 
     252      NAMELIST/namsbc_alb/ nn_ice_alb, rn_albice  
    203253      !!---------------------------------------------------------------------- 
    204254      ! 
     
    219269         WRITE(numout,*) '~~~~~~~' 
    220270         WRITE(numout,*) '   Namelist namsbc_alb : albedo ' 
    221          WRITE(numout,*) '      correction for snow and ice albedo                  rn_cloud  = ', rn_cloud 
    222          WRITE(numout,*) '      albedo of melting ice in the arctic and antarctic   rn_albice = ', rn_albice 
    223          WRITE(numout,*) '      coefficients for linear                             rn_alphd  = ', rn_alphd 
    224          WRITE(numout,*) '      interpolation used to compute albedo                rn_alphdi = ', rn_alphdi 
    225          WRITE(numout,*) '      between two extremes values (Pyane, 1972)           rn_alphc  = ', rn_alphc 
     271         WRITE(numout,*) '      choose the albedo parameterization                  nn_ice_alb = ', nn_ice_alb 
     272         WRITE(numout,*) '      albedo of bare puddled ice                          rn_albice  = ', rn_albice 
    226273      ENDIF 
    227274      ! 
Note: See TracChangeset for help on using the changeset viewer.