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 9937 for NEMO/branches/2018/dev_r9838_ENHANCE04_MLF/src/OCE – NEMO

Ignore:
Timestamp:
2018-07-12T17:55:41+02:00 (6 years ago)
Author:
gm
Message:

#1911 (ENHANCE-04): step I.2 (end): clean sea ice related physical constant in dev_r9838_ENHANCE04_MLF

Location:
NEMO/branches/2018/dev_r9838_ENHANCE04_MLF/src/OCE
Files:
7 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2018/dev_r9838_ENHANCE04_MLF/src/OCE/BDY/bdyice.F90

    r9923 r9937  
    124124 
    125125            ! Then, a) transfer the snow excess into the ice (different from icethd_dh) 
    126             zdh = MAX( 0._wp, ( rhosn * h_s(ji,jj,jl) + ( rhoic - rho0 ) * h_i(ji,jj,jl) ) * r1_rho0 ) 
     126            zdh = MAX( 0._wp, ( rhos * h_s(ji,jj,jl) + ( rhoi - rho0 ) * h_i(ji,jj,jl) ) * r1_rho0 ) 
    127127            ! Or, b) transfer all the snow into ice (if incoming ice is likely to melt as it comes into a warmer environment) 
    128             !zdh = MAX( 0._wp, h_s(ji,jj,jl) * rhosn / rhoic ) 
     128            !zdh = MAX( 0._wp, h_s(ji,jj,jl) * rhos / rhoi ) 
    129129 
    130130            ! recompute h_i, h_s 
    131131            h_i(ji,jj,jl) = MIN( hi_max(jl), h_i(ji,jj,jl) + zdh ) 
    132             h_s(ji,jj,jl) = MAX( 0._wp, h_s(ji,jj,jl) - zdh * rhoic / rhosn )  
    133  
    134          ENDDO 
     132            h_s(ji,jj,jl) = MAX( 0._wp, h_s(ji,jj,jl) - zdh * rhoi / rhos )  
     133 
     134         END DO 
    135135         CALL lbc_bdy_lnk( a_i(:,:,jl), 'T', 1., ib_bdy ) 
    136136         CALL lbc_bdy_lnk( h_i(:,:,jl), 'T', 1., ib_bdy ) 
    137137         CALL lbc_bdy_lnk( h_s(:,:,jl), 'T', 1., ib_bdy ) 
    138       ENDDO 
     138      END DO 
    139139      ! retrieve at_i 
    140140      at_i(:,:) = 0._wp 
     
    212212            DO jk = 1, nlay_s 
    213213               ! Snow energy of melting 
    214                e_s(ji,jj,jk,jl) = rswitch * rhosn * ( cpic * ( rt0 - t_s(ji,jj,jk,jl) ) + lfus ) 
     214               e_s(ji,jj,jk,jl) = rswitch * rhos * ( rcpi * ( rt0 - t_s(ji,jj,jk,jl) ) + rLfus ) 
    215215               ! Multiply by volume, so that heat content in J/m2 
    216216               e_s(ji,jj,jk,jl) = e_s(ji,jj,jk,jl) * v_s(ji,jj,jl) * r1_nlay_s 
     
    219219               ztmelts          = - tmut * sz_i(ji,jj,jk,jl) + rt0 !Melting temperature in K                   
    220220               ! heat content per unit volume 
    221                e_i(ji,jj,jk,jl) = rswitch * rhoic * & 
    222                   (   cpic    * ( ztmelts - t_i(ji,jj,jk,jl) ) & 
    223                   +   lfus    * ( 1.0 - (ztmelts-rt0) / MIN((t_i(ji,jj,jk,jl)-rt0),-epsi20) ) & 
    224                   - rcp      * ( ztmelts - rt0 ) ) 
     221               e_i(ji,jj,jk,jl) = rswitch * rhoi * & 
     222                  (   rcpi  * ( ztmelts - t_i(ji,jj,jk,jl) ) & 
     223                  +   rLfus * ( 1.0 - (ztmelts-rt0) / MIN((t_i(ji,jj,jk,jl)-rt0),-epsi20) ) & 
     224                  -   rcp   * ( ztmelts - rt0 ) ) 
    225225               ! Mutliply by ice volume, and divide by number of layers to get heat content in J/m2 
    226226               e_i(ji,jj,jk,jl) = e_i(ji,jj,jk,jl) * a_i(ji,jj,jl) * h_i(ji,jj,jl) * r1_nlay_i 
  • NEMO/branches/2018/dev_r9838_ENHANCE04_MLF/src/OCE/DOM/phycst.F90

    r9923 r9937  
    3434   REAL(wp), PUBLIC ::   rhhmm =  60._wp        !: number of minutes in one hour 
    3535   REAL(wp), PUBLIC ::   rmmss =  60._wp        !: number of seconds in one minute 
    36    REAL(wp), PUBLIC ::   omega                  !: earth rotation parameter           [s-1] 
    37    REAL(wp), PUBLIC ::   ra    = 6371229._wp    !: earth radius                       [m] 
    38    REAL(wp), PUBLIC ::   grav  = 9.80665_wp     !: gravity                            [m/s2] 
     36   REAL(wp), PUBLIC ::   omega                  !: earth rotation parameter               [s-1] 
     37   REAL(wp), PUBLIC ::   ra    = 6371229._wp    !: earth radius                           [m] 
     38   REAL(wp), PUBLIC ::   grav  = 9.80665_wp     !: gravity                                [m/s2] 
    3939    
    40    REAL(wp), PUBLIC ::   rtt      = 273.16_wp        !: triple point of temperature   [Kelvin] 
    41    REAL(wp), PUBLIC ::   rt0      = 273.15_wp        !: freezing point of fresh water [Kelvin] 
    42    REAL(wp), PUBLIC ::   rt0_snow = 273.15_wp        !: melting point of snow         [Kelvin] 
    43 #if defined key_si3 
    44    REAL(wp), PUBLIC ::   rt0_ice  = 273.15_wp        !: melting point of ice          [Kelvin] 
    45 #else 
    46    REAL(wp), PUBLIC ::   rt0_ice  = 273.05_wp        !: melting point of ice          [Kelvin] 
    47 #endif 
    48    REAL(wp), PUBLIC ::   rho0                        !: volumic mass of reference     [kg/m3] 
    49    REAL(wp), PUBLIC ::   r1_rho0                     !: = 1. / rho0                   [m3/kg] 
    50    REAL(wp), PUBLIC ::   rcp                         !: ocean specific heat           [J/Kelvin] 
    51    REAL(wp), PUBLIC ::   r1_rcp                      !: = 1. / rcp                    [Kelvin/J] 
     40   REAL(wp), PUBLIC ::   rt0      = 273.15_wp        !: freezing point of fresh water     [Kelvin] 
     41   REAL(wp), PUBLIC ::   rho0                        !: volumic mass of reference         [kg/m3] 
     42   REAL(wp), PUBLIC ::   r1_rho0                     !: = 1. / rho0                       [m3/kg] 
     43   REAL(wp), PUBLIC ::   rcp                         !: ocean specific heat               [J/Kelvin] 
     44   REAL(wp), PUBLIC ::   r1_rcp                      !: = 1. / rcp                        [Kelvin/J] 
    5245   REAL(wp), PUBLIC ::   rho0_rcp                    !: = rho0 * rcp  
    5346   REAL(wp), PUBLIC ::   r1_rho0_rcp                 !: = 1. / ( rho0 * rcp ) 
    5447 
    55    REAL(wp), PUBLIC ::   rhosn    =  330._wp         !: volumic mass of snow          [kg/m3] 
    56    REAL(wp), PUBLIC ::   rhofw    = 1000._wp         !: volumic mass of freshwater in melt ponds [kg/m3] 
     48   REAL(wp), PUBLIC ::   rhoi     =  917._wp         !: sea ice density                   [kg/m3] 
     49   REAL(wp), PUBLIC ::   rhos     =  330._wp         !: snow    density                   [kg/m3] 
     50   REAL(wp), PUBLIC ::   rhow     = 1000._wp         !: water   density (in melt ponds)   [kg/m3] 
     51   REAL(wp), PUBLIC ::   rcnd_i   =    2.034396_wp   !: thermal conductivity of fresh ice [W/m/K] 
     52   REAL(wp), PUBLIC ::   rcpi     = 2067.0_wp        !: specific heat of fresh ice        [J/kg/K] 
     53   REAL(wp), PUBLIC ::   rLsub    =    2.834e+6_wp   !: pure ice latent heat of sublimation   [J/kg] 
     54   REAL(wp), PUBLIC ::   rLfus    =    0.334e+6_wp   !: latent heat of fusion of fresh ice    [J/kg] 
     55   REAL(wp), PUBLIC ::   tmut     =    0.054_wp      !: decrease of seawater meltpoint with salinity 
    5756   REAL(wp), PUBLIC ::   emic     =    0.97_wp       !: emissivity of snow or ice 
    58    REAL(wp), PUBLIC ::   sice     =    6.0_wp        !: salinity of ice               [psu] 
    59    REAL(wp), PUBLIC ::   soce     =   34.7_wp        !: salinity of sea               [psu] 
    60    REAL(wp), PUBLIC ::   cevap    =    2.5e+6_wp     !: latent heat of evaporation (water) 
    61    REAL(wp), PUBLIC ::   srgamma  =    0.9_wp        !: correction factor for solar radiation (Oberhuber, 1974) 
     57   REAL(wp), PUBLIC ::   sice     =    6.0_wp        !: salinity of ice                   [psu] 
     58   REAL(wp), PUBLIC ::   soce     =   34.7_wp        !: salinity of sea                   [psu] 
    6259   REAL(wp), PUBLIC ::   vkarmn   =    0.4_wp        !: von Karman constant 
    6360   REAL(wp), PUBLIC ::   stefan   =    5.67e-8_wp    !: Stefan-Boltzmann constant  
    6461 
    65 #if defined key_si3 || defined key_cice 
    66    REAL(wp), PUBLIC ::   rhoic    =  917._wp         !: volumic mass of sea ice                               [kg/m3] 
    67    REAL(wp), PUBLIC ::   rcdic    =    2.034396_wp   !: thermal conductivity of fresh ice                     [W/m/K] 
    68    REAL(wp), PUBLIC ::   cpic     = 2067.0_wp        !: specific heat of fresh ice                            [J/kg/K] 
    69    REAL(wp), PUBLIC ::   lsub     =    2.834e+6_wp   !: pure ice latent heat of sublimation                   [J/kg] 
    70    REAL(wp), PUBLIC ::   lfus     =    0.334e+6_wp   !: latent heat of fusion of fresh ice                    [J/kg] 
    71    REAL(wp), PUBLIC ::   tmut     =    0.054_wp      !: decrease of seawater meltpoint with salinity 
    72    REAL(wp), PUBLIC ::   xlsn                        !: = lfus*rhosn (volumetric latent heat fusion of snow)  [J/m3] 
    73 #else 
    74    REAL(wp), PUBLIC ::   rhoic    =  900._wp         !: volumic mass of sea ice                               [kg/m3] 
    75    REAL(wp), PUBLIC ::   rcdic    =    2.034396_wp   !: conductivity of the ice                               [W/m/K] 
    76    REAL(wp), PUBLIC ::   rcpic    =    1.8837e+6_wp  !: volumetric specific heat for ice                      [J/m3/K] 
    77    REAL(wp), PUBLIC ::   cpic                        !: = rcpic / rhoic  (specific heat for ice)              [J/Kg/K] 
    78    REAL(wp), PUBLIC ::   rcdsn    =    0.22_wp       !: conductivity of the snow                              [W/m/K] 
    79    REAL(wp), PUBLIC ::   rcpsn    =    6.9069e+5_wp  !: volumetric specific heat for snow                     [J/m3/K] 
    80    REAL(wp), PUBLIC ::   xlsn     =  110.121e+6_wp   !: volumetric latent heat fusion of snow                 [J/m3] 
    81    REAL(wp), PUBLIC ::   lfus                        !: = xlsn / rhosn   (latent heat of fusion of fresh ice) [J/Kg] 
    82    REAL(wp), PUBLIC ::   xlic     =  300.33e+6_wp    !: volumetric latent heat fusion of ice                  [J/m3] 
    83    REAL(wp), PUBLIC ::   xsn      =    2.8e+6_wp     !: volumetric latent heat of sublimation of snow         [J/m3] 
    84 #endif 
    85 #if defined key_cice 
    86    REAL(wp), PUBLIC ::   rcdsn    =    0.31_wp       !: thermal conductivity of snow                          [W/m/K] 
    87 #endif 
    88 #if defined key_si3 
    89    REAL(wp), PUBLIC ::   r1_rhoic                    !: 1 / rhoic 
    90    REAL(wp), PUBLIC ::   r1_rhosn                    !: 1 / rhosn 
    91    REAL(wp), PUBLIC ::   r1_cpic                     !: 1 / cpic 
    92 #endif 
     62   REAL(wp), PUBLIC ::   r1_rhoi                     !: 1 / rhoi 
     63   REAL(wp), PUBLIC ::   r1_rhos                     !: 1 / rhos 
     64   REAL(wp), PUBLIC ::   r1_rhow                     !: 1 / rhow 
     65   REAL(wp), PUBLIC ::   r1_cpi                      !: 1 / rcpi 
     66   REAL(wp), PUBLIC ::   r1_Lsub                     !: 1 / rLsub 
     67   REAL(wp), PUBLIC ::   r1_Lfus                     !: 1 / rLfus 
     68 
    9369   !!---------------------------------------------------------------------- 
    9470   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    10581      !! ** Purpose :   set and print the constants 
    10682      !!---------------------------------------------------------------------- 
    107  
     83      ! 
    10884      IF(lwp) WRITE(numout,*) 
    10985      IF(lwp) WRITE(numout,*) 'phy_cst : initialization of ocean parameters and constants' 
    11086      IF(lwp) WRITE(numout,*) '~~~~~~~' 
    11187 
    112       ! Define & print constants 
    113       ! ------------------------ 
    114       IF(lwp) WRITE(numout,*) 
    115       IF(lwp) WRITE(numout,*) '   Constants' 
    116  
    117       IF(lwp) WRITE(numout,*) 
    118       IF(lwp) WRITE(numout,*) '      mathematical constant                 rpi = ', rpi 
     88      !                 !==  Define derived constant  ==! 
    11989 
    12090      rsiyea = 365.25_wp * rday * 2._wp * rpi / 6.283076_wp 
     
    12595      omega  = 2._wp * rpi / rsiday  
    12696#endif 
    127       IF(lwp) WRITE(numout,*) 
    128       IF(lwp) WRITE(numout,*) '      day                                rday   = ', rday,   ' s' 
    129       IF(lwp) WRITE(numout,*) '      sideral year                       rsiyea = ', rsiyea, ' s' 
    130       IF(lwp) WRITE(numout,*) '      sideral day                        rsiday = ', rsiday, ' s' 
    131       IF(lwp) WRITE(numout,*) '      omega                              omega  = ', omega,  ' s^-1' 
    132       IF(lwp) WRITE(numout,*) 
    133       IF(lwp) WRITE(numout,*) '      nb of months per year               raamo = ', raamo, ' months' 
    134       IF(lwp) WRITE(numout,*) '      nb of hours per day                 rjjhh = ', rjjhh, ' hours' 
    135       IF(lwp) WRITE(numout,*) '      nb of minutes per hour              rhhmm = ', rhhmm, ' mn' 
    136       IF(lwp) WRITE(numout,*) '      nb of seconds per minute            rmmss = ', rmmss, ' s' 
    137       IF(lwp) WRITE(numout,*) 
    138       IF(lwp) WRITE(numout,*) '      earth radius                         ra   = ', ra, ' m' 
    139       IF(lwp) WRITE(numout,*) '      gravity                              grav = ', grav , ' m/s^2' 
    140       IF(lwp) WRITE(numout,*) 
    141       IF(lwp) WRITE(numout,*) '      triple point of temperature      rtt      = ', rtt     , ' K' 
    142       IF(lwp) WRITE(numout,*) '      freezing point of water          rt0      = ', rt0     , ' K' 
    143       IF(lwp) WRITE(numout,*) '      melting point of snow            rt0_snow = ', rt0_snow, ' K' 
    144       IF(lwp) WRITE(numout,*) '      melting point of ice             rt0_ice  = ', rt0_ice , ' K' 
    145       IF(lwp) WRITE(numout,*) 
    146       IF(lwp) WRITE(numout,*) '   reference density and heat capacity now defined in eosbn2.f90' 
    147                
    148 #if defined key_si3 || defined key_cice 
    149       xlsn = lfus * rhosn        ! volumetric latent heat fusion of snow [J/m3] 
    150 #else 
    151       cpic = rcpic / rhoic       ! specific heat for ice   [J/Kg/K] 
    152       lfus = xlsn / rhosn        ! latent heat of fusion of fresh ice 
    153 #endif 
    154 #if defined key_si3 
    155       r1_rhoic = 1._wp / rhoic 
    156       r1_rhosn = 1._wp / rhosn 
    157       r1_cpic  = 1._wp / cpic 
    158 #endif 
    159       IF(lwp) THEN 
     97 
     98      r1_rhoi = 1._wp / rhoi 
     99      r1_rhos = 1._wp / rhos 
     100      r1_cpi  = 1._wp / rcpi 
     101      r1_Lsub = 1._wp / rLsub 
     102      r1_Lfus = 1._wp / rLfus 
     103 
     104      IF(lwp) THEN      !==  print constants  ==! 
    160105         WRITE(numout,*) 
    161 #if defined key_cice 
    162          WRITE(numout,*) '      thermal conductivity of the snow          = ', rcdsn   , ' J/s/m/K' 
    163 #endif 
    164          WRITE(numout,*) '      thermal conductivity of pure ice          = ', rcdic   , ' J/s/m/K' 
    165          WRITE(numout,*) '      fresh ice specific heat                   = ', cpic    , ' J/kg/K' 
    166          WRITE(numout,*) '      latent heat of fusion of fresh ice / snow = ', lfus    , ' J/kg' 
    167 #if defined key_si3 || defined key_cice 
    168          WRITE(numout,*) '      latent heat of subl.  of fresh ice / snow = ', lsub    , ' J/kg' 
    169 #else 
    170          WRITE(numout,*) '      density times specific heat for snow      = ', rcpsn   , ' J/m^3/K'  
    171          WRITE(numout,*) '      density times specific heat for ice       = ', rcpic   , ' J/m^3/K' 
    172          WRITE(numout,*) '      volumetric latent heat fusion of sea ice  = ', xlic    , ' J/m'  
    173          WRITE(numout,*) '      latent heat of sublimation of snow        = ', xsn     , ' J/kg'  
    174 #endif 
    175          WRITE(numout,*) '      volumetric latent heat fusion of snow     = ', xlsn    , ' J/m^3'  
    176          WRITE(numout,*) '      density of sea ice                        = ', rhoic   , ' kg/m^3' 
    177          WRITE(numout,*) '      density of snow                           = ', rhosn   , ' kg/m^3' 
    178          WRITE(numout,*) '      density of freshwater (in melt ponds)     = ', rhofw   , ' kg/m^3' 
    179          WRITE(numout,*) '      emissivity of snow or ice                 = ', emic   
    180          WRITE(numout,*) '      salinity of ice                           = ', sice    , ' psu' 
    181          WRITE(numout,*) '      salinity of sea                           = ', soce    , ' psu' 
    182          WRITE(numout,*) '      latent heat of evaporation (water)        = ', cevap   , ' J/m^3'  
    183          WRITE(numout,*) '      correction factor for solar radiation     = ', srgamma  
    184          WRITE(numout,*) '      von Karman constant                       = ', vkarmn  
    185          WRITE(numout,*) '      Stefan-Boltzmann constant                 = ', stefan  , ' J/s/m^2/K^4' 
     106         WRITE(numout,*) '   Constants' 
    186107         WRITE(numout,*) 
    187          WRITE(numout,*) '      conversion: degre ==> radian          rad = ', rad 
     108         WRITE(numout,*) '      mathematical constant              rpi    = ', rpi 
     109         WRITE(numout,*) '      conversion: degre ==> radian       rad    = ', rad 
    188110         WRITE(numout,*) 
    189          WRITE(numout,*) '      smallest real computer value       rsmall = ', rsmall 
     111         WRITE(numout,*) '      day in seconds                     rday   = ', rday  , ' s' 
     112         WRITE(numout,*) '      sideral year                       rsiyea = ', rsiyea, ' s' 
     113         WRITE(numout,*) '      sideral day                        rsiday = ', rsiday, ' s' 
     114         WRITE(numout,*) '      omega = 2 pi / rsiday              omega  = ', omega , ' s^-1' 
     115         WRITE(numout,*) '      earth radius                       ra     = ', ra    , ' m' 
     116         WRITE(numout,*) '      gravity                            grav   = ', grav  , ' m/s^2' 
     117         WRITE(numout,*) 
     118         WRITE(numout,*) '      nb of months per year              raamo  = ', raamo, ' months' 
     119         WRITE(numout,*) '      nb of hours per day                rjjhh  = ', rjjhh, ' hours' 
     120         WRITE(numout,*) '      nb of minutes per hour             rhhmm  = ', rhhmm, ' mn' 
     121         WRITE(numout,*) '      nb of seconds per minute           rmmss  = ', rmmss, ' s' 
     122         WRITE(numout,*) 
     123         WRITE(numout,*) '   reference ocean density and heat capacity now defined in eosbn2.f90' 
     124         WRITE(numout,*) 
     125         WRITE(numout,*) '      freezing point of freshwater                rt0    = ', rt0   , ' K' 
     126         WRITE(numout,*) '      sea ice density                             rhoi   = ', rhoi  , ' kg/m^3' 
     127         WRITE(numout,*) '      snow    density                             rhos   = ', rhos  , ' kg/m^3' 
     128         WRITE(numout,*) '      freshwater density (in melt ponds)          rhow   = ', rhow  , ' kg/m^3' 
     129         WRITE(numout,*) '      thermal conductivity of pure ice            rcnd_i = ', rcnd_i, ' J/s/m/K' 
     130         WRITE(numout,*) '      fresh ice specific heat                     rcpi   = ', rcpi  , ' J/kg/K' 
     131         WRITE(numout,*) '      latent heat of fusion of fresh ice / snow   rLfus  = ', rLfus , ' J/kg' 
     132         WRITE(numout,*) '      latent heat of subl.  of fresh ice / snow   rLsub  = ', rLsub , ' J/kg' 
     133         WRITE(numout,*) '      emissivity of snow or ice                   emic   = ', emic   
     134         WRITE(numout,*) '      salinity of ice                             sice   = ', sice  , ' psu' 
     135         WRITE(numout,*) '      salinity of sea                             soce   = ', soce  , ' psu' 
     136         WRITE(numout,*) '      von Karman constant                         vkarmn = ', vkarmn  
     137         WRITE(numout,*) '      Stefan-Boltzmann constant                   stefan = ', stefan, ' J/s/m^2/K^4' 
     138         WRITE(numout,*) 
     139         WRITE(numout,*) 
     140         WRITE(numout,*) '      smallest real computer value                rsmall = ', rsmall 
    190141      ENDIF 
    191  
     142      ! 
    192143   END SUBROUTINE phy_cst 
    193144 
  • NEMO/branches/2018/dev_r9838_ENHANCE04_MLF/src/OCE/SBC/sbcblk.F90

    r9923 r9937  
    504504      ENDIF 
    505505 
    506       zqla(:,:) = L_vap(zst(:,:)) * zevap(:,:)     ! Latent Heat flux 
     506      zqla(:,:) = L_vap( zst(:,:) ) * zevap(:,:)     ! Latent Heat flux 
    507507 
    508508 
     
    526526      ! 
    527527      qns(:,:) = zqlw(:,:) - zqsb(:,:) - zqla(:,:)                                &   ! Downward Non Solar 
    528          &     - sf(jp_snow)%fnow(:,:,1) * rn_pfac * lfus                         &   ! remove latent melting heat for solid precip 
     528         &     - sf(jp_snow)%fnow(:,:,1) * rn_pfac * rLfus                        &   ! remove latent melting heat for solid precip 
    529529         &     - zevap(:,:) * pst(:,:) * rcp                                      &   ! remove evap heat content at SST 
    530530         &     + ( sf(jp_prec)%fnow(:,:,1) - sf(jp_snow)%fnow(:,:,1) ) * rn_pfac  &   ! add liquid precip heat content at Tair 
    531531         &     * ( sf(jp_tair)%fnow(:,:,1) - rt0 ) * rcp                          & 
    532532         &     + sf(jp_snow)%fnow(:,:,1) * rn_pfac                                &   ! add solid  precip heat content at min(Tair,Tsnow) 
    533          &     * ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic 
     533         &     * ( MIN( sf(jp_tair)%fnow(:,:,1), rt0 ) - rt0 ) * rcpi 
    534534      qns(:,:) = qns(:,:) * tmask(:,:,1) 
    535535      ! 
     
    643643      !! ** Purpose : Compute the moist adiabatic lapse-rate. 
    644644      !!     => http://glossary.ametsoc.org/wiki/Moist-adiabatic_lapse_rate 
    645       !!     => http://www.geog.ucsb.edu/~joel/g266_s10/lecture_notes/chapt03/oh10_3_01/oh10_3_01.html 
    646645      !! 
    647646      !! ** Author: L. Brodeau, june 2016 / AeroBulk (https://sourceforge.net/p/aerobulk) 
     
    652651      ! 
    653652      INTEGER  ::   ji, jj         ! dummy loop indices 
    654       REAL(wp) :: zrv, ziRT        ! local scalar 
     653      REAL(wp) ::   zrv, ziRT      ! local scalar 
     654      REAL(wp) ::   zLv = 2.5e+6_wp   ! latent heat of vaporisation  
    655655      !!---------------------------------------------------------------------------------- 
    656656      ! 
     
    659659            zrv = pqa(ji,jj) / (1. - pqa(ji,jj)) 
    660660            ziRT = 1. / (R_dry*ptak(ji,jj))    ! 1/RT 
    661             gamma_moist(ji,jj) = grav * ( 1. + cevap*zrv*ziRT ) / ( Cp_dry + cevap*cevap*zrv*reps0*ziRT/ptak(ji,jj) ) 
     661            gamma_moist(ji,jj) = grav * ( 1. + zLv*zrv*ziRT ) / ( Cp_dry + zLv*zLv*zrv*reps0*ziRT/ptak(ji,jj) ) 
    662662         END DO 
    663663      END DO 
     
    792792      REAL(wp) ::   zst3                     ! local variable 
    793793      REAL(wp) ::   zcoef_dqlw, zcoef_dqla   !   -      - 
    794       REAL(wp) ::   zztmp, z1_lsub           !   -      - 
     794      REAL(wp) ::   zztmp                    !   -      - 
    795795      REAL(wp) ::   zfr1, zfr2               ! local variables 
    796796      REAL(wp), DIMENSION(jpi,jpj,jpl) ::   z1_st         ! inverse of surface temperature 
     
    868868 
    869869      ! --- evaporation --- ! 
    870       z1_lsub = 1._wp / Lsub 
    871       evap_ice (:,:,:) = rn_efac * qla_ice (:,:,:) * z1_lsub    ! sublimation 
    872       devap_ice(:,:,:) = rn_efac * dqla_ice(:,:,:) * z1_lsub    ! d(sublimation)/dT 
     870      evap_ice (:,:,:) = rn_efac * qla_ice (:,:,:) * r1_Lsub    ! sublimation 
     871      devap_ice(:,:,:) = rn_efac * dqla_ice(:,:,:) * r1_Lsub    ! d(sublimation)/dT 
    873872      zevap    (:,:)   = rn_efac * ( emp(:,:) + tprecip(:,:) )  ! evaporation over ocean 
    874873 
     
    884883         &          + ( tprecip(:,:) - sprecip(:,:) ) * ( sf(jp_tair)%fnow(:,:,1) - rt0 ) * rcp  & ! liquid precip at Tair 
    885884         &          +   sprecip(:,:) * ( 1._wp - zsnw ) *                                        & ! solid precip at min(Tair,Tsnow) 
    886          &              ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus ) 
     885         &              ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0 ) - rt0 ) * rcpi * tmask(:,:,1) - rLfus ) 
    887886      qemp_ice(:,:) =   sprecip(:,:) * zsnw *                                                    & ! solid precip (only) 
    888          &              ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus ) 
     887         &              ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0 ) - rt0 ) * rcpi * tmask(:,:,1) - rLfus ) 
    889888 
    890889      ! --- total solar and non solar fluxes --- ! 
     
    894893 
    895894      ! --- heat content of precip over ice in J/m3 (to be used in 1D-thermo) --- ! 
    896       qprec_ice(:,:) = rhosn * ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0_snow ) - rt0 ) * cpic * tmask(:,:,1) - lfus ) 
     895      qprec_ice(:,:) = rhos * ( ( MIN( sf(jp_tair)%fnow(:,:,1), rt0 ) - rt0 ) * rcpi * tmask(:,:,1) - rLfus ) 
    897896 
    898897      ! --- heat content of evap over ice in W/m2 (to be used in 1D-thermo) --- 
    899898      DO jl = 1, jpl 
    900          qevap_ice(:,:,jl) = 0._wp ! should be -evap_ice(:,:,jl)*( ( Tice - rt0 ) * cpic * tmask(:,:,1) ) 
     899         qevap_ice(:,:,jl) = 0._wp ! should be -evap_ice(:,:,jl)*( ( Tice - rt0 ) * rcpi * tmask(:,:,1) ) 
    901900         !                         ! But we do not have Tice => consider it at 0degC => evap=0  
    902901      END DO 
     
    971970      CASE ( 1 , 2 ) 
    972971         ! 
    973          zfac  = 1._wp /  ( rn_cnd_s + rcdic ) 
     972         zfac  = 1._wp /  ( rn_cnd_s + rcnd_i ) 
    974973         zfac2 = EXP(1._wp) * 0.5_wp * zepsilon 
    975974         zfac3 = 2._wp / zepsilon 
     
    978977            DO jj = 1 , jpj 
    979978               DO ji = 1, jpi 
    980                   zhe = ( rn_cnd_s * phi(ji,jj,jl) + rcdic * phs(ji,jj,jl) ) * zfac                             ! Effective thickness 
    981                   IF( zhe >=  zfac2 )   zgfac(ji,jj,jl) = MIN( 2._wp, 0.5_wp * ( 1._wp + LOG( zhe * zfac3 ) ) ) ! Enhanced conduction factor 
     979                  zhe = ( rn_cnd_s * phi(ji,jj,jl) + rcnd_i * phs(ji,jj,jl) ) * zfac                             ! Effective thickness 
     980                  IF( zhe >=  zfac2 )   zgfac(ji,jj,jl) = MIN( 2._wp, 0.5_wp * ( 1._wp + LOG( zhe * zfac3 ) ) )  ! Enhanced conduction factor 
    982981               END DO 
    983982            END DO 
     
    990989      ! -------------------------------------------------------------! 
    991990      ! 
    992       zfac = rcdic * rn_cnd_s 
     991      zfac = rcnd_i * rn_cnd_s 
    993992      ! 
    994993      DO jl = 1, jpl 
    995994         DO jj = 1 , jpj 
    996995            DO ji = 1, jpi 
    997                !                     
    998                zkeff_h = zfac * zgfac(ji,jj,jl) / &                                    ! Effective conductivity of the snow-ice system divided by thickness 
    999                   &      ( rcdic * phs(ji,jj,jl) + rn_cnd_s * MAX( 0.01, phi(ji,jj,jl) ) ) 
     996               !                                                                       ! Effective conductivity of the snow-ice system divided by thickness 
     997               zkeff_h = zfac * zgfac(ji,jj,jl) / ( rcnd_i * phs(ji,jj,jl) + rn_cnd_s * MAX( 0.01, phi(ji,jj,jl) ) ) 
    1000998               ztsu    = ptsu(ji,jj,jl)                                                ! Store current iteration temperature 
    1001999               ztsu0   = ptsu(ji,jj,jl)                                                ! Store initial surface temperature 
  • NEMO/branches/2018/dev_r9838_ENHANCE04_MLF/src/OCE/SBC/sbccpl.F90

    r9923 r9937  
    14181418            zqns(:,:) =  zqns(:,:) - zemp(:,:) * sst_m(:,:) * rcp         ! remove heat content due to mass flux (assumed to be at SST) 
    14191419            IF( srcv(jpr_snow  )%laction ) THEN 
    1420                zqns(:,:) = zqns(:,:) - frcv(jpr_snow)%z3(:,:,1) * lfus    ! energy for melting solid precipitation over the free ocean 
     1420               zqns(:,:) = zqns(:,:) - frcv(jpr_snow)%z3(:,:,1) * rLfus   ! energy for melting solid precipitation over the free ocean 
    14211421            ENDIF 
    14221422         ENDIF 
    14231423         ! 
    1424          IF( srcv(jpr_icb)%laction )  zqns(:,:) = zqns(:,:) - frcv(jpr_icb)%z3(:,:,1) * lfus ! remove heat content associated to iceberg melting 
     1424         IF( srcv(jpr_icb)%laction )  zqns(:,:) = zqns(:,:) - frcv(jpr_icb)%z3(:,:,1) * rLfus ! remove heat content associated to iceberg melting 
    14251425         ! 
    14261426         IF( ln_mixcpl ) THEN   ;   qns(:,:) = qns(:,:) * xcplmask(:,:,0) + zqns(:,:) * zmsk(:,:) 
     
    18111811      !                                      
    18121812      ! --- calving (removed from qns_tot) --- ! 
    1813       IF( srcv(jpr_cal)%laction )   zqns_tot(:,:) = zqns_tot(:,:) - frcv(jpr_cal)%z3(:,:,1) * lfus  ! remove latent heat of calving 
    1814                                                                                                     ! we suppose it melts at 0deg, though it should be temp. of surrounding ocean 
     1813      IF( srcv(jpr_cal)%laction )   zqns_tot(:,:) = zqns_tot(:,:) - frcv(jpr_cal)%z3(:,:,1) * rLfus  ! remove latent heat of calving 
     1814                                                                                                     ! we suppose it melts at 0deg, though it should be temp. of surrounding ocean 
    18151815      ! --- iceberg (removed from qns_tot) --- ! 
    1816       IF( srcv(jpr_icb)%laction )   zqns_tot(:,:) = zqns_tot(:,:) - frcv(jpr_icb)%z3(:,:,1) * lfus  ! remove latent heat of iceberg melting 
     1816      IF( srcv(jpr_icb)%laction )   zqns_tot(:,:) = zqns_tot(:,:) - frcv(jpr_icb)%z3(:,:,1) * rLfus  ! remove latent heat of iceberg melting 
    18171817 
    18181818#if defined key_si3       
     
    18231823 
    18241824      ! Heat content per unit mass of snow (J/kg) 
    1825       WHERE( SUM( a_i, dim=3 ) > 1.e-10 )   ;   zcptsnw(:,:) = cpic * SUM( (tn_ice - rt0) * a_i, dim=3 ) / SUM( a_i, dim=3 ) 
     1825      WHERE( SUM( a_i, dim=3 ) > 1.e-10 )   ;   zcptsnw(:,:) = rcpi * SUM( (tn_ice - rt0) * a_i, dim=3 ) / SUM( a_i, dim=3 ) 
    18261826      ELSEWHERE                             ;   zcptsnw(:,:) = zcptn(:,:) 
    1827       ENDWHERE 
     1827      END WHERE 
    18281828      ! Heat content per unit mass of rain (J/kg) 
    18291829      zcptrain(:,:) = rcp * ( SUM( (tn_ice(:,:,:) - rt0) * a_i(:,:,:), dim=3 ) + sst_m(:,:) * ziceld(:,:) )  
    18301830 
    18311831      ! --- enthalpy of snow precip over ice in J/m3 (to be used in 1D-thermo) --- ! 
    1832       zqprec_ice(:,:) = rhosn * ( zcptsnw(:,:) - lfus ) 
     1832      zqprec_ice(:,:) = rhos * ( zcptsnw(:,:) - rLfus ) 
    18331833 
    18341834      ! --- heat content of evap over ice in W/m2 (to be used in 1D-thermo) --- ! 
    18351835      DO jl = 1, jpl 
    1836          zqevap_ice(:,:,jl) = 0._wp ! should be -evap * ( ( Tice - rt0 ) * cpic ) but atm. does not take it into account 
     1836         zqevap_ice(:,:,jl) = 0._wp ! should be -evap * ( ( Tice - rt0 ) * rcpi ) but atm. does not take it into account 
    18371837      END DO 
    18381838 
     
    18401840      zqemp_oce(:,:) = -  zevap_oce(:,:)                                      *   zcptn   (:,:)   &        ! evap 
    18411841         &             + ( ztprecip(:,:) - zsprecip(:,:) )                    *   zcptrain(:,:)   &        ! liquid precip 
    1842          &             +   zsprecip(:,:)                   * ( 1._wp - zsnw ) * ( zcptsnw (:,:) - lfus )   ! solid precip over ocean + snow melting 
    1843       zqemp_ice(:,:) =     zsprecip(:,:)                   * zsnw             * ( zcptsnw (:,:) - lfus )   ! solid precip over ice (qevap_ice=0 since atm. does not take it into account) 
     1842         &             +   zsprecip(:,:)                   * ( 1._wp - zsnw ) * ( zcptsnw (:,:) - rLfus )   ! solid precip over ocean + snow melting 
     1843      zqemp_ice(:,:) =     zsprecip(:,:)                   * zsnw             * ( zcptsnw (:,:) - rLfus )   ! solid precip over ice (qevap_ice=0 since atm. does not take it into account) 
    18441844!!    zqemp_ice(:,:) = -   frcv(jpr_ievp)%z3(:,:,1)        * picefr(:,:)      *   zcptsnw (:,:)   &        ! ice evap 
    1845 !!       &             +   zsprecip(:,:)                   * zsnw             * zqprec_ice(:,:) * r1_rhosn ! solid precip over ice 
     1845!!       &             +   zsprecip(:,:)                   * zsnw             * zqprec_ice(:,:) * r1_rhos  ! solid precip over ice 
    18461846       
    18471847      ! --- total non solar flux (including evap/precip) --- ! 
     
    18751875      ! clem: this formulation is certainly wrong... but better than it was... 
    18761876      zqns_tot(:,:) = zqns_tot(:,:)                            &          ! zqns_tot update over free ocean with: 
    1877          &          - (  ziceld(:,:) * zsprecip(:,:) * lfus )  &          ! remove the latent heat flux of solid precip. melting 
     1877         &          - (  ziceld(:,:) * zsprecip(:,:) * rLfus )  &         ! remove the latent heat flux of solid precip. melting 
    18781878         &          - (  zemp_tot(:,:)                         &          ! remove the heat content of mass flux (assumed to be at SST) 
    18791879         &             - zemp_ice(:,:) ) * zcptn(:,:)  
     
    18921892#endif 
    18931893      ! outputs 
    1894       IF ( srcv(jpr_cal)%laction       ) CALL iom_put('hflx_cal_cea'    , - frcv(jpr_cal)%z3(:,:,1) * lfus )                       ! latent heat from calving 
    1895       IF ( srcv(jpr_icb)%laction       ) CALL iom_put('hflx_icb_cea'    , - frcv(jpr_icb)%z3(:,:,1) * lfus )                       ! latent heat from icebergs melting 
     1894      IF ( srcv(jpr_cal)%laction       ) CALL iom_put('hflx_cal_cea'    , - frcv(jpr_cal)%z3(:,:,1) * rLfus )                      ! latent heat from calving 
     1895      IF ( srcv(jpr_icb)%laction       ) CALL iom_put('hflx_icb_cea'    , - frcv(jpr_icb)%z3(:,:,1) * rLfus )                      ! latent heat from icebergs melting 
    18961896      IF ( iom_use('hflx_rain_cea')    ) CALL iom_put('hflx_rain_cea'   , ( tprecip(:,:) - sprecip(:,:) ) * zcptrain(:,:) )        ! heat flux from rain (cell average) 
    18971897      IF ( iom_use('hflx_evap_cea')    ) CALL iom_put('hflx_evap_cea'   , ( frcv(jpr_tevp)%z3(:,:,1) - frcv(jpr_ievp)%z3(:,:,1) & 
    18981898           &                                                              * picefr(:,:) ) * zcptn(:,:) * tmask(:,:,1) )            ! heat flux from evap (cell average) 
    1899       IF ( iom_use('hflx_snow_cea')    ) CALL iom_put('hflx_snow_cea'   , sprecip(:,:) * ( zcptsnw(:,:) - Lfus )   )               ! heat flux from snow (cell average) 
    1900       IF ( iom_use('hflx_snow_ao_cea') ) CALL iom_put('hflx_snow_ao_cea', sprecip(:,:) * ( zcptsnw(:,:) - Lfus ) & 
     1899      IF ( iom_use('hflx_snow_cea')    ) CALL iom_put('hflx_snow_cea'   , sprecip(:,:) * ( zcptsnw(:,:) - rLfus )   )               ! heat flux from snow (cell average) 
     1900      IF ( iom_use('hflx_snow_ao_cea') ) CALL iom_put('hflx_snow_ao_cea', sprecip(:,:) * ( zcptsnw(:,:) - rLfus ) & 
    19011901           &                                                              * ( 1._wp - zsnw(:,:) )                  )               ! heat flux from snow (over ocean) 
    1902       IF ( iom_use('hflx_snow_ai_cea') ) CALL iom_put('hflx_snow_ai_cea', sprecip(:,:) * ( zcptsnw(:,:) - Lfus ) & 
     1902      IF ( iom_use('hflx_snow_ai_cea') ) CALL iom_put('hflx_snow_ai_cea', sprecip(:,:) * ( zcptsnw(:,:) - rLfus ) & 
    19031903           &                                                              *           zsnw(:,:)                    )               ! heat flux from snow (over ice) 
    19041904      ! note: hflx for runoff and iceshelf are done in sbcrnf and sbcisf resp. 
  • NEMO/branches/2018/dev_r9838_ENHANCE04_MLF/src/OCE/SBC/sbcice_cice.F90

    r9923 r9937  
    1313   USE dom_oce         ! ocean space and time domain 
    1414   USE domvvl 
    15    USE phycst   , ONLY : rcp, rho0, r1_rho0, rhosn, rhoic 
     15   USE phycst   , ONLY : rcp, rho0, r1_rho0, rhos, rhoi 
    1616   USE in_out_manager  ! I/O manager 
    1717   USE iom, ONLY : iom_put,iom_use              ! I/O manager library !!Joakim edit 
     
    222222      CALL cice2nemo(vsno(:,:,:),ztmp1,'T', 1. ) 
    223223      CALL cice2nemo(vice(:,:,:),ztmp2,'T', 1. ) 
    224       snwice_mass  (:,:) = ( rhosn * ztmp1(:,:) + rhoic * ztmp2(:,:)  ) 
     224      snwice_mass  (:,:) = ( rhos * ztmp1(:,:) + rhoi * ztmp2(:,:)  ) 
    225225      snwice_mass_b(:,:) = snwice_mass(:,:) 
    226226 
     
    644644      CALL cice2nemo(vsno(:,:,:),ztmp1,'T', 1. ) 
    645645      CALL cice2nemo(vice(:,:,:),ztmp2,'T', 1. ) 
    646       snwice_mass  (:,:) = ( rhosn * ztmp1(:,:) + rhoic * ztmp2(:,:)  ) 
     646      snwice_mass  (:,:) = ( rhos * ztmp1(:,:) + rhoi * ztmp2(:,:)  ) 
    647647      snwice_mass_b(:,:) = snwice_mass(:,:) 
    648648      snwice_fmass (:,:) = ( snwice_mass(:,:) - snwice_mass_b(:,:) ) / dt 
  • NEMO/branches/2018/dev_r9838_ENHANCE04_MLF/src/OCE/SBC/sbcisf.F90

    r9923 r9937  
    5252   LOGICAL, PUBLIC ::   l_isfcpl = .false.       !: isf recieved from oasis 
    5353 
    54    REAL(wp), PUBLIC, SAVE ::   rcpi     = 2000.0_wp     !: specific heat of ice shelf             [J/kg/K] 
     54   REAL(wp)        , SAVE ::   rcp_isf  = 2000.0_wp     !: specific heat of ice shelf             [J/kg/K] 
    5555   REAL(wp), PUBLIC, SAVE ::   rkappa   = 1.54e-6_wp    !: heat diffusivity through the ice-shelf [m2/s] 
    5656   REAL(wp), PUBLIC, SAVE ::   rho_isf  = 920.0_wp      !: volumic mass of ice shelf              [kg/m3] 
    5757   REAL(wp), PUBLIC, SAVE ::   tsurf    = -20.0_wp      !: air temperature on top of ice shelf    [C] 
    58    REAL(wp), PUBLIC, SAVE ::   rlfusisf = 0.334e6_wp    !: latent heat of fusion of ice shelf     [J/kg] 
    5958 
    6059!: Variable used in fldread to read the forcing file (nn_isf == 4 .OR. nn_isf == 3) 
     
    114113            ! compute fwf and heat flux 
    115114            IF( .NOT.l_isfcpl ) THEN    ;   CALL sbc_isf_cav (kt) 
    116             ELSE                        ;   qisf(:,:)  = fwfisf(:,:) * rlfusisf  ! heat        flux 
     115            ELSE                        ;   qisf(:,:)  = fwfisf(:,:) * rLfus    ! heat flux 
    117116            ENDIF 
    118117            ! 
     
    127126               fwfisf(:,:) = - sf_rnfisf(1)%fnow(:,:,1)         ! fresh water flux from the isf (fwfisf <0 mean melting)  
    128127            ENDIF 
    129             qisf(:,:)   = fwfisf(:,:) * rlfusisf             ! heat flux 
     128            qisf(:,:)   = fwfisf(:,:) * rLfus                   ! heat flux 
    130129            stbl(:,:)   = soce 
    131130            ! 
     
    137136               fwfisf(:,:) = -sf_fwfisf(1)%fnow(:,:,1)            ! fwf 
    138137            ENDIF 
    139             qisf(:,:)   = fwfisf(:,:) * rlfusisf               ! heat flux 
     138            qisf(:,:)   = fwfisf(:,:) * rLfus                     ! heat flux 
    140139            stbl(:,:)   = soce 
    141140            ! 
     
    308307      qisf    (:,:)    = 0._wp   ;   fwfisf  (:,:) = 0._wp 
    309308      risf_tsc(:,:,:)  = 0._wp   ;   fwfisf_b(:,:) = 0._wp 
    310       ! 
    311       ! define isf tbl tickness, top and bottom indice 
    312       SELECT CASE ( nn_isf ) 
     309 
     310      SELECT CASE ( nn_isf )      ! define isf tbl tickness, top and bottom indice 
     311      ! 
    313312      CASE ( 1 )  
    314313         IF(lwp) WRITE(numout,*) 
     
    455454                           & * r1_e1e2t(ji,jj) * tmask(ji,jj,jk) 
    456455              
    457                fwfisf(ji,jj) = qisf(ji,jj) / rlfusisf          !fresh water flux kg/(m2s)                   
     456               fwfisf(ji,jj) = qisf(ji,jj) * r1_Lfus                        ! fresh water flux kg/(m2s)                   
    458457               fwfisf(ji,jj) = fwfisf(ji,jj) * ( soce / stbl(ji,jj) ) 
    459458               !add to salinity trend 
     
    527526               DO ji = 1, jpi 
    528527                  zhtflx(ji,jj) =   zgammat(ji,jj)*rcp*rho0*(ttbl(ji,jj)-zfrz(ji,jj)) 
    529                   zfwflx(ji,jj) = - zhtflx(ji,jj)/rlfusisf 
     528                  zfwflx(ji,jj) = - zhtflx(ji,jj) * r1_Lfus 
    530529               END DO 
    531530            END DO 
     
    545544                  ! compute coeficient to solve the 2nd order equation 
    546545                  zeps1 = rcp*rho0*zgammat(ji,jj) 
    547                   zeps2 = rlfusisf*rho0*zgammas(ji,jj) 
    548                   zeps3 = rho_isf*rcpi*rkappa/MAX(risfdep(ji,jj),zeps) 
     546                  zeps2 = rLfus*rho0*zgammas(ji,jj) 
     547                  zeps3 = rho_isf*rcp_isf*rkappa/MAX(risfdep(ji,jj),zeps) 
    549548                  zeps4 = zlamb2+zlamb3*risfdep(ji,jj) 
    550549                  zeps6 = zeps4-ttbl(ji,jj) 
  • NEMO/branches/2018/dev_r9838_ENHANCE04_MLF/src/OCE/SBC/sbcrnf.F90

    r9923 r9937  
    128128            END WHERE 
    129129            WHERE( sf_t_rnf(1)%fnow(:,:,1) == -222._wp )             ! where fwf comes from melting of ice shelves or iceberg 
    130                rnf_tsc(:,:,jp_tem) = ztfrz(:,:) * rnf(:,:) * r1_rho0 - rnf(:,:) * rlfusisf * r1_rho0_rcp 
     130               rnf_tsc(:,:,jp_tem) = ztfrz(:,:) * rnf(:,:) * r1_rho0 - rnf(:,:) * rLfus * r1_rho0_rcp 
    131131            END WHERE 
    132132         ELSE                                                        ! use SST as runoffs temperature 
Note: See TracChangeset for help on using the changeset viewer.