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.
limthd_dh.F90 in trunk/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: trunk/NEMOGCM/NEMO/LIM_SRC_3/limthd_dh.F90 @ 5202

Last change on this file since 5202 was 5202, checked in by clem, 9 years ago

LIM3: important bug fix to avoid crashes. See ticket #1508

  • Property svn:keywords set to Id
File size: 35.8 KB
RevLine 
[825]1MODULE limthd_dh
[1572]2   !!======================================================================
3   !!                       ***  MODULE limthd_dh ***
4   !!  LIM-3 :   thermodynamic growth and decay of the ice
5   !!======================================================================
6   !! History :  LIM  ! 2003-05 (M. Vancoppenolle) Original code in 1D
7   !!                 ! 2005-06 (M. Vancoppenolle) 3D version
[4688]8   !!            3.2  ! 2009-07 (M. Vancoppenolle, Y. Aksenov, G. Madec) bug correction in wfx_snw & wfx_ice
[3625]9   !!            3.4  ! 2011-02 (G. Madec) dynamical allocation
10   !!            3.5  ! 2012-10 (G. Madec & co) salt flux + bug fixes
[1572]11   !!----------------------------------------------------------------------
[825]12#if defined key_lim3
[834]13   !!----------------------------------------------------------------------
14   !!   'key_lim3'                                      LIM3 sea-ice model
15   !!----------------------------------------------------------------------
[3625]16   !!   lim_thd_dh    : vertical accr./abl. and lateral ablation of sea ice
[825]17   !!----------------------------------------------------------------------
[3625]18   USE par_oce        ! ocean parameters
19   USE phycst         ! physical constants (OCE directory)
20   USE sbc_oce        ! Surface boundary condition: ocean fields
21   USE ice            ! LIM variables
22   USE thd_ice        ! LIM thermodynamics
23   USE in_out_manager ! I/O manager
24   USE lib_mpp        ! MPP library
25   USE wrk_nemo       ! work arrays
26   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
[4688]27   
[825]28   IMPLICIT NONE
29   PRIVATE
30
[1572]31   PUBLIC   lim_thd_dh   ! called by lim_thd
[825]32
33   !!----------------------------------------------------------------------
[4161]34   !! NEMO/LIM3 4.0 , UCL - NEMO Consortium (2010)
[1156]35   !! $Id$
[2715]36   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[825]37   !!----------------------------------------------------------------------
38CONTAINS
39
[4688]40   SUBROUTINE lim_thd_dh( kideb, kiut )
[921]41      !!------------------------------------------------------------------
42      !!                ***  ROUTINE lim_thd_dh  ***
43      !!
[1572]44      !! ** Purpose :   determines variations of ice and snow thicknesses.
[921]45      !!
[1572]46      !! ** Method  :   Ice/Snow surface melting arises from imbalance in surface fluxes
47      !!              Bottom accretion/ablation arises from flux budget
48      !!              Snow thickness can increase by precipitation and decrease by sublimation
49      !!              If snow load excesses Archmiede limit, snow-ice is formed by
50      !!              the flooding of sea-water in the snow
[921]51      !!
[1572]52      !!                 1) Compute available flux of heat for surface ablation
53      !!                 2) Compute snow and sea ice enthalpies
54      !!                 3) Surface ablation and sublimation
55      !!                 4) Bottom accretion/ablation
56      !!                 5) Case of Total ablation
57      !!                 6) Snow ice formation
[921]58      !!
[1572]59      !! References : Bitz and Lipscomb, 1999, J. Geophys. Res.
60      !!              Fichefet T. and M. Maqueda 1997, J. Geophys. Res., 102(C6), 12609-12646   
61      !!              Vancoppenolle, Fichefet and Bitz, 2005, Geophys. Res. Let.
62      !!              Vancoppenolle et al.,2009, Ocean Modelling
[921]63      !!------------------------------------------------------------------
[1572]64      INTEGER , INTENT(in) ::   kideb, kiut   ! Start/End point on which the  the computation is applied
65      !!
66      INTEGER  ::   ji , jk        ! dummy loop indices
[4161]67      INTEGER  ::   ii, ij         ! 2D corresponding indices to ji
[1572]68      INTEGER  ::   iter
[825]69
[4688]70      REAL(wp) ::   ztmelts             ! local scalar
[5123]71      REAL(wp) ::   zfdum       
[1572]72      REAL(wp) ::   zfracs       ! fractionation coefficient for bottom salt entrapment
73      REAL(wp) ::   zcoeff       ! dummy argument for snowfall partitioning over ice and leads
[4688]74      REAL(wp) ::   zs_snic  ! snow-ice salinity
[1572]75      REAL(wp) ::   zswi1        ! switch for computation of bottom salinity
76      REAL(wp) ::   zswi12       ! switch for computation of bottom salinity
77      REAL(wp) ::   zswi2        ! switch for computation of bottom salinity
78      REAL(wp) ::   zgrr         ! bottom growth rate
[4688]79      REAL(wp) ::   zt_i_new     ! bottom formation temperature
80
81      REAL(wp) ::   zQm          ! enthalpy exchanged with the ocean (J/m2), >0 towards the ocean
82      REAL(wp) ::   zEi          ! specific enthalpy of sea ice (J/kg)
83      REAL(wp) ::   zEw          ! specific enthalpy of exchanged water (J/kg)
84      REAL(wp) ::   zdE          ! specific enthalpy difference (J/kg)
85      REAL(wp) ::   zfmdt        ! exchange mass flux x time step (J/m2), >0 towards the ocean
86      REAL(wp) ::   zsstK        ! SST in Kelvin
87
88      REAL(wp), POINTER, DIMENSION(:) ::   zqprec      ! energy of fallen snow                       (J.m-3)
89      REAL(wp), POINTER, DIMENSION(:) ::   zq_su       ! heat for surface ablation                   (J.m-2)
90      REAL(wp), POINTER, DIMENSION(:) ::   zq_bo       ! heat for bottom ablation                    (J.m-2)
91      REAL(wp), POINTER, DIMENSION(:) ::   zq_rema     ! remaining heat at the end of the routine    (J.m-2)
[5123]92      REAL(wp), POINTER, DIMENSION(:) ::   zf_tt       ! Heat budget to determine melting or freezing(W.m-2)
[3294]93
[3625]94      REAL(wp), POINTER, DIMENSION(:) ::   zdh_s_mel   ! snow melt
95      REAL(wp), POINTER, DIMENSION(:) ::   zdh_s_pre   ! snow precipitation
96      REAL(wp), POINTER, DIMENSION(:) ::   zdh_s_sub   ! snow sublimation
[3294]97
98      REAL(wp), POINTER, DIMENSION(:,:) ::   zdeltah
[4688]99      REAL(wp), POINTER, DIMENSION(:,:) ::   zh_i      ! ice layer thickness
[5202]100      INTEGER , POINTER, DIMENSION(:,:) ::   icount    ! number of layers vanished by melting
[3294]101
[4688]102      REAL(wp), POINTER, DIMENSION(:) ::   zqh_i       ! total ice heat content  (J.m-2)
103      REAL(wp), POINTER, DIMENSION(:) ::   zqh_s       ! total snow heat content (J.m-2)
104      REAL(wp), POINTER, DIMENSION(:) ::   zq_s        ! total snow enthalpy     (J.m-3)
[3294]105
[5123]106      REAL(wp) :: zswitch_sal
[4161]107
[3294]108      ! Heat conservation
[4688]109      INTEGER  ::   num_iter_max
110
[1572]111      !!------------------------------------------------------------------
[825]112
[5123]113      ! Discriminate between varying salinity (nn_icesal=2) and prescribed cases (other values)
114      SELECT CASE( nn_icesal )                       ! varying salinity or not
[4688]115         CASE( 1, 3, 4 ) ;   zswitch_sal = 0       ! prescribed salinity profile
116         CASE( 2 )       ;   zswitch_sal = 1       ! varying salinity profile
117      END SELECT
[825]118
[5167]119      CALL wrk_alloc( jpij, zqprec, zq_su, zq_bo, zf_tt, zq_rema )
[4688]120      CALL wrk_alloc( jpij, zdh_s_mel, zdh_s_pre, zdh_s_sub, zqh_i, zqh_s, zq_s )
[5202]121      CALL wrk_alloc( jpij, nlay_i, zdeltah, zh_i )
122      CALL wrk_alloc( jpij, nlay_i, icount )
[4161]123     
[4688]124      dh_i_surf  (:) = 0._wp ; dh_i_bott  (:) = 0._wp ; dh_snowice(:) = 0._wp
125      dsm_i_se_1d(:) = 0._wp ; dsm_i_si_1d(:) = 0._wp   
126 
127      zqprec (:) = 0._wp ; zq_su  (:) = 0._wp ; zq_bo  (:) = 0._wp ; zf_tt  (:) = 0._wp
[5123]128      zq_rema(:) = 0._wp
[2715]129
[4688]130      zdh_s_pre(:) = 0._wp
131      zdh_s_mel(:) = 0._wp
132      zdh_s_sub(:) = 0._wp
133      zqh_s    (:) = 0._wp     
134      zqh_i    (:) = 0._wp   
[4161]135
[4688]136      zh_i      (:,:) = 0._wp       
137      zdeltah   (:,:) = 0._wp       
[5202]138      icount    (:,:) = 0
[4688]139
[5167]140      ! Initialize enthalpy at nlay_i+1
141      DO ji = kideb, kiut
142         q_i_1d(ji,nlay_i+1) = 0._wp
143      END DO
144
[4688]145      ! initialize layer thicknesses and enthalpies
146      h_i_old (:,0:nlay_i+1) = 0._wp
147      qh_i_old(:,0:nlay_i+1) = 0._wp
148      DO jk = 1, nlay_i
149         DO ji = kideb, kiut
[5123]150            h_i_old (ji,jk) = ht_i_1d(ji) * r1_nlay_i
[4872]151            qh_i_old(ji,jk) = q_i_1d(ji,jk) * h_i_old(ji,jk)
[4688]152         ENDDO
153      ENDDO
[921]154      !
155      !------------------------------------------------------------------------------!
[4688]156      !  1) Calculate available heat for surface and bottom ablation                 !
[921]157      !------------------------------------------------------------------------------!
158      !
[2715]159      DO ji = kideb, kiut
[4990]160         zfdum      = qns_ice_1d(ji) + ( 1._wp - i0(ji) ) * qsr_ice_1d(ji) - fc_su(ji) 
161         zf_tt(ji)  = fc_bo_i(ji) + fhtur_1d(ji) + fhld_1d(ji) 
[4688]162
[5167]163         zq_su (ji) = MAX( 0._wp, zfdum     * rdt_ice ) * MAX( 0._wp , SIGN( 1._wp, t_su_1d(ji) - rt0 ) )
[4688]164         zq_bo (ji) = MAX( 0._wp, zf_tt(ji) * rdt_ice )
165      END DO
166
[921]167      !
168      !------------------------------------------------------------------------------!
[4688]169      ! If snow temperature is above freezing point, then snow melts
170      ! (should not happen but sometimes it does)
[921]171      !------------------------------------------------------------------------------!
[4688]172      DO ji = kideb, kiut
[5123]173         IF( t_s_1d(ji,1) > rt0 ) THEN !!! Internal melting
[4688]174            ! Contribution to heat flux to the ocean [W.m-2], < 0 
[4872]175            hfx_res_1d(ji) = hfx_res_1d(ji) + q_s_1d(ji,1) * ht_s_1d(ji) * a_i_1d(ji) * r1_rdtice
[4688]176            ! Contribution to mass flux
[4872]177            wfx_snw_1d(ji) = wfx_snw_1d(ji) + rhosn * ht_s_1d(ji) * a_i_1d(ji) * r1_rdtice
[4688]178            ! updates
[4872]179            ht_s_1d(ji)   = 0._wp
180            q_s_1d (ji,1) = 0._wp
[5123]181            t_s_1d (ji,1) = rt0
[4688]182         END IF
183      END DO
184
185      !------------------------------------------------------------!
186      !  2) Computing layer thicknesses and enthalpies.            !
187      !------------------------------------------------------------!
[921]188      !
[825]189      DO jk = 1, nlay_s
[2715]190         DO ji = kideb, kiut
[5167]191            zqh_s(ji) =  zqh_s(ji) + q_s_1d(ji,jk) * ht_s_1d(ji) * r1_nlay_s
[825]192         END DO
193      END DO
[2715]194      !
[825]195      DO jk = 1, nlay_i
[2715]196         DO ji = kideb, kiut
[5123]197            zh_i(ji,jk) = ht_i_1d(ji) * r1_nlay_i
[4872]198            zqh_i(ji)   = zqh_i(ji) + q_i_1d(ji,jk) * zh_i(ji,jk)
[825]199         END DO
200      END DO
[921]201      !
202      !------------------------------------------------------------------------------|
203      !  3) Surface ablation and sublimation                                         |
204      !------------------------------------------------------------------------------|
205      !
[834]206      !-------------------------
207      ! 3.1 Snow precips / melt
208      !-------------------------
[825]209      ! Snow accumulation in one thermodynamic time step
210      ! snowfall is partitionned between leads and ice
211      ! if snow fall was uniform, a fraction (1-at_i) would fall into leads
212      ! but because of the winds, more snow falls on leads than on sea ice
213      ! and a greater fraction (1-at_i)^beta of the total mass of snow
[834]214      ! (beta < 1) falls in leads.
[825]215      ! In reality, beta depends on wind speed,
216      ! and should decrease with increasing wind speed but here, it is
[834]217      ! considered as a constant. an average value is 0.66
[825]218      ! Martin Vancoppenolle, December 2006
219
[5167]220      zdeltah(:,:) = 0._wp
[825]221      DO ji = kideb, kiut
[4688]222         !-----------
223         ! Snow fall
224         !-----------
225         ! thickness change
[5123]226         zcoeff = ( 1._wp - ( 1._wp - at_i_1d(ji) )**rn_betas ) / at_i_1d(ji) 
227         zdh_s_pre(ji) = zcoeff * sprecip_1d(ji) * rdt_ice * r1_rhosn
[4688]228         ! enthalpy of the precip (>0, J.m-3) (tatm_ice is now in K)
[5123]229         zqprec   (ji) = rhosn * ( cpic * ( rt0 - MIN( tatm_ice_1d(ji), rt0_snow) ) + lfus )   
[4688]230         IF( sprecip_1d(ji) == 0._wp ) zqprec(ji) = 0._wp
231         ! heat flux from snow precip (>0, W.m-2)
[4872]232         hfx_spr_1d(ji) = hfx_spr_1d(ji) + zdh_s_pre(ji) * a_i_1d(ji) * zqprec(ji) * r1_rdtice
[4688]233         ! mass flux, <0
[4872]234         wfx_spr_1d(ji) = wfx_spr_1d(ji) - rhosn * a_i_1d(ji) * zdh_s_pre(ji) * r1_rdtice
[825]235
[4688]236         !---------------------
237         ! Melt of falling snow
238         !---------------------
239         ! thickness change
[5134]240         rswitch        = MAX( 0._wp , SIGN( 1._wp , zqprec(ji) - epsi20 ) )
[5167]241         zdeltah (ji,1) = - rswitch * zq_su(ji) / MAX( zqprec(ji) , epsi20 )
242         zdeltah (ji,1) = MAX( - zdh_s_pre(ji), zdeltah(ji,1) ) ! bound melting
[4688]243         ! heat used to melt snow (W.m-2, >0)
[5167]244         hfx_snw_1d(ji) = hfx_snw_1d(ji) - zdeltah(ji,1) * a_i_1d(ji) * zqprec(ji) * r1_rdtice
[4688]245         ! snow melting only = water into the ocean (then without snow precip), >0
[5167]246         wfx_snw_1d(ji) = wfx_snw_1d(ji) - rhosn * a_i_1d(ji) * zdeltah(ji,1) * r1_rdtice   
247         ! updates available heat + precipitations after melting
248         zq_su     (ji) = MAX( 0._wp , zq_su (ji) + zdeltah(ji,1) * zqprec(ji) )     
249         zdh_s_pre (ji) = zdh_s_pre(ji) + zdeltah(ji,1)
[4688]250
[5167]251         ! update thickness
252         ht_s_1d(ji) = MAX( 0._wp , ht_s_1d(ji) + zdh_s_pre(ji) )
[825]253      END DO
254
[5167]255      ! If heat still available (zq_su > 0), then melt more snow
256      zdeltah(:,:) = 0._wp
[825]257      DO jk = 1, nlay_s
258         DO ji = kideb, kiut
[4688]259            ! thickness change
[4990]260            rswitch          = 1._wp - MAX( 0._wp, SIGN( 1._wp, - ht_s_1d(ji) ) ) 
[5134]261            rswitch          = rswitch * ( MAX( 0._wp, SIGN( 1._wp, q_s_1d(ji,jk) - epsi20 ) ) ) 
[4990]262            zdeltah  (ji,jk) = - rswitch * zq_su(ji) / MAX( q_s_1d(ji,jk), epsi20 )
[5167]263            zdeltah  (ji,jk) = MAX( zdeltah(ji,jk) , - ht_s_1d(ji) ) ! bound melting
[4688]264            zdh_s_mel(ji)    = zdh_s_mel(ji) + zdeltah(ji,jk)   
265            ! heat used to melt snow(W.m-2, >0)
[4872]266            hfx_snw_1d(ji)   = hfx_snw_1d(ji) - zdeltah(ji,jk) * a_i_1d(ji) * q_s_1d(ji,jk) * r1_rdtice 
[4688]267            ! snow melting only = water into the ocean (then without snow precip)
[4872]268            wfx_snw_1d(ji)   = wfx_snw_1d(ji) - rhosn * a_i_1d(ji) * zdeltah(ji,jk) * r1_rdtice
[4688]269            ! updates available heat + thickness
[5123]270            zq_su (ji)  = MAX( 0._wp , zq_su (ji) + zdeltah(ji,jk) * q_s_1d(ji,jk) )
[4872]271            ht_s_1d(ji) = MAX( 0._wp , ht_s_1d(ji) + zdeltah(ji,jk) )
[1572]272         END DO
273      END DO
[825]274
[4688]275      !----------------------
276      ! 3.2 Snow sublimation
277      !----------------------
278      ! qla_ice is always >=0 (upwards), heat goes to the atmosphere, therefore snow sublimates
[5167]279      ! clem comment: not counted in mass/heat exchange in limsbc since this is an exchange with atm. (not ocean)
[4688]280      ! clem comment: ice should also sublimate
[5167]281      zdeltah(:,:) = 0._wp
[4688]282      IF( lk_cpl ) THEN
283         ! coupled mode: sublimation already included in emp_ice (to do in limsbc_ice)
284         zdh_s_sub(:)      =  0._wp 
285      ELSE
286         ! forced  mode: snow thickness change due to sublimation
[921]287         DO ji = kideb, kiut
[5167]288            zdh_s_sub(ji)  =  MAX( - ht_s_1d(ji) , - qla_ice_1d(ji) / ( rhosn * lsub ) * rdt_ice )
[4688]289            ! Heat flux by sublimation [W.m-2], < 0
290            !      sublimate first snow that had fallen, then pre-existing snow
[5167]291            zdeltah(ji,1)  = MAX( zdh_s_sub(ji), - zdh_s_pre(ji) )
292            hfx_sub_1d(ji) = hfx_sub_1d(ji) + ( zdeltah(ji,1) * zqprec(ji) + ( zdh_s_sub(ji) - zdeltah(ji,1) ) * q_s_1d(ji,1)  &
293               &                              ) * a_i_1d(ji) * r1_rdtice
[4688]294            ! Mass flux by sublimation
[4872]295            wfx_sub_1d(ji) =  wfx_sub_1d(ji) - rhosn * a_i_1d(ji) * zdh_s_sub(ji) * r1_rdtice
[4688]296            ! new snow thickness
[5167]297            ht_s_1d(ji)    =  MAX( 0._wp , ht_s_1d(ji) + zdh_s_sub(ji) )
298            ! update precipitations after sublimation and correct sublimation
299            zdh_s_pre(ji) = zdh_s_pre(ji) + zdeltah(ji,1)
300            zdh_s_sub(ji) = zdh_s_sub(ji) - zdeltah(ji,1)
[1572]301         END DO
302      ENDIF
[825]303
[4688]304      ! --- Update snow diags --- !
[825]305      DO ji = kideb, kiut
[5167]306         dh_s_tot(ji) = zdh_s_mel(ji) + zdh_s_pre(ji) + zdh_s_sub(ji)
307      END DO
[825]308
[4688]309      !-------------------------------------------
310      ! 3.3 Update temperature, energy
311      !-------------------------------------------
312      ! new temp and enthalpy of the snow (remaining snow precip + remaining pre-existing snow)
313      zq_s(:) = 0._wp 
[825]314      DO jk = 1, nlay_s
315         DO ji = kideb,kiut
[5134]316            rswitch       = MAX(  0._wp , SIGN( 1._wp, ht_s_1d(ji) - epsi20 )  )
317            q_s_1d(ji,jk) = rswitch / MAX( ht_s_1d(ji), epsi20 ) *                          &
[5167]318              &            ( (   zdh_s_pre(ji)             ) * zqprec(ji) +  &
319              &              (   ht_s_1d(ji) - zdh_s_pre(ji) ) * rhosn * ( cpic * ( rt0 - t_s_1d(ji,jk) ) + lfus ) )
[4872]320            zq_s(ji)     =  zq_s(ji) + q_s_1d(ji,jk)
[825]321         END DO
322      END DO
323
[4688]324      !--------------------------
325      ! 3.4 Surface ice ablation
326      !--------------------------
327      zdeltah(:,:) = 0._wp ! important
328      DO jk = 1, nlay_i
[5167]329         DO ji = kideb, kiut
[5202]330            ztmelts           = - tmut * s_i_1d(ji,jk) + rt0          ! Melting point of layer k [K]
331           
332            IF( t_i_1d(ji,jk) >= ztmelts ) THEN !!! Internal melting
[4688]333
[5202]334               zEi            = - q_i_1d(ji,jk) * r1_rhoic            ! Specific enthalpy of layer k [J/kg, <0]       
335               zdE            = 0._wp                                 ! Specific enthalpy difference   (J/kg, <0)
336                                                                      ! set up at 0 since no energy is needed to melt water...(it is already melted)
337               zdeltah(ji,jk) = MIN( 0._wp , - zh_i(ji,jk) )          ! internal melting occurs when the internal temperature is above freezing     
338                                                                      ! this should normally not happen, but sometimes, heat diffusion leads to this
339               zfmdt          = - zdeltah(ji,jk) * rhoic              ! Mass flux x time step > 0
340                         
341               dh_i_surf(ji)  = dh_i_surf(ji) + zdeltah(ji,jk)        ! Cumulate surface melt
342               
343               zfmdt          = - rhoic * zdeltah(ji,jk)              ! Recompute mass flux [kg/m2, >0]
[4688]344
[5202]345               ! Contribution to heat flux to the ocean [W.m-2], <0 (ice enthalpy zEi is "sent" to the ocean)
346               hfx_res_1d(ji) = hfx_res_1d(ji) + zfmdt * a_i_1d(ji) * zEi * r1_rdtice
347               
348               ! Contribution to salt flux (clem: using sm_i_1d and not s_i_1d(jk) is ok)
349               sfx_res_1d(ji) = sfx_res_1d(ji) - rhoic * a_i_1d(ji) * zdeltah(ji,jk) * sm_i_1d(ji) * r1_rdtice
350               
351               ! Contribution to mass flux
352               wfx_res_1d(ji) =  wfx_res_1d(ji) - rhoic * a_i_1d(ji) * zdeltah(ji,jk) * r1_rdtice
[4688]353
[5202]354            ELSE                                !!! Surface melting
355               
356               zEi            = - q_i_1d(ji,jk) * r1_rhoic            ! Specific enthalpy of layer k [J/kg, <0]
357               zEw            =    rcp * ( ztmelts - rt0 )            ! Specific enthalpy of resulting meltwater [J/kg, <0]
358               zdE            =    zEi - zEw                          ! Specific enthalpy difference < 0
359               
360               zfmdt          = - zq_su(ji) / zdE                     ! Mass flux to the ocean [kg/m2, >0]
361               
362               zdeltah(ji,jk) = - zfmdt * r1_rhoic                    ! Melt of layer jk [m, <0]
363               
364               zdeltah(ji,jk) = MIN( 0._wp , MAX( zdeltah(ji,jk) , - zh_i(ji,jk) ) )    ! Melt of layer jk cannot exceed the layer thickness [m, <0]
365               
366               zq_su(ji)      = MAX( 0._wp , zq_su(ji) - zdeltah(ji,jk) * rhoic * zdE ) ! update available heat
367               
368               dh_i_surf(ji)  = dh_i_surf(ji) + zdeltah(ji,jk)        ! Cumulate surface melt
369               
370               zfmdt          = - rhoic * zdeltah(ji,jk)              ! Recompute mass flux [kg/m2, >0]
371               
372               zQm            = zfmdt * zEw                           ! Energy of the melt water sent to the ocean [J/m2, <0]
373               
374               ! Contribution to salt flux (clem: using sm_i_1d and not s_i_1d(jk) is ok)
375               sfx_sum_1d(ji) = sfx_sum_1d(ji) - rhoic * a_i_1d(ji) * zdeltah(ji,jk) * sm_i_1d(ji) * r1_rdtice
376               
377               ! Contribution to heat flux [W.m-2], < 0
378               hfx_thd_1d(ji) = hfx_thd_1d(ji) + zfmdt * a_i_1d(ji) * zEw * r1_rdtice
379               
380               ! Total heat flux used in this process [W.m-2], > 0 
381               hfx_sum_1d(ji) = hfx_sum_1d(ji) - zfmdt * a_i_1d(ji) * zdE * r1_rdtice
382               
383               ! Contribution to mass flux
384               wfx_sum_1d(ji) =  wfx_sum_1d(ji) - rhoic * a_i_1d(ji) * zdeltah(ji,jk) * r1_rdtice
385               
[5167]386            END IF
[4688]387            ! record which layers have disappeared (for bottom melting)
388            !    => icount=0 : no layer has vanished
389            !    => icount=5 : 5 layers have vanished
[5202]390            rswitch       = MAX( 0._wp , SIGN( 1._wp , - ( zh_i(ji,jk) + zdeltah(ji,jk) ) ) ) 
391            icount(ji,jk) = NINT( rswitch )
392            zh_i(ji,jk)   = MAX( 0._wp , zh_i(ji,jk) + zdeltah(ji,jk) )
[4688]393
394            ! update heat content (J.m-2) and layer thickness
[4872]395            qh_i_old(ji,jk) = qh_i_old(ji,jk) + zdeltah(ji,jk) * q_i_1d(ji,jk)
[4688]396            h_i_old (ji,jk) = h_i_old (ji,jk) + zdeltah(ji,jk)
[825]397         END DO
[921]398      END DO
[4688]399      ! update ice thickness
400      DO ji = kideb, kiut
[4872]401         ht_i_1d(ji) =  MAX( 0._wp , ht_i_1d(ji) + dh_i_surf(ji) )
[4688]402      END DO
[825]403
[921]404      !
405      !------------------------------------------------------------------------------!
406      ! 4) Basal growth / melt                                                       !
407      !------------------------------------------------------------------------------!
408      !
[4688]409      !------------------
410      ! 4.1 Basal growth
411      !------------------
412      ! Basal growth is driven by heat imbalance at the ice-ocean interface,
413      ! between the inner conductive flux  (fc_bo_i), from the open water heat flux
414      ! (fhld) and the turbulent ocean flux (fhtur).
415      ! fc_bo_i is positive downwards. fhtur and fhld are positive to the ice
[825]416
[4688]417      ! If salinity varies in time, an iterative procedure is required, because
418      ! the involved quantities are inter-dependent.
419      ! Basal growth (dh_i_bott) depends upon new ice specific enthalpy (zEi),
420      ! which depends on forming ice salinity (s_i_new), which depends on dh/dt (dh_i_bott)
421      ! -> need for an iterative procedure, which converges quickly
422
[5167]423      num_iter_max = 1
424      IF( nn_icesal == 2 ) num_iter_max = 5
[825]425
[4688]426      ! Iterative procedure
427      DO iter = 1, num_iter_max
428         DO ji = kideb, kiut
429            IF(  zf_tt(ji) < 0._wp  ) THEN
[825]430
[4688]431               ! New bottom ice salinity (Cox & Weeks, JGR88 )
432               !--- zswi1  if dh/dt < 2.0e-8
433               !--- zswi12 if 2.0e-8 < dh/dt < 3.6e-7
434               !--- zswi2  if dh/dt > 3.6e-7
435               zgrr               = MIN( 1.0e-3, MAX ( dh_i_bott(ji) * r1_rdtice , epsi10 ) )
436               zswi2              = MAX( 0._wp , SIGN( 1._wp , zgrr - 3.6e-7 ) )
437               zswi12             = MAX( 0._wp , SIGN( 1._wp , zgrr - 2.0e-8 ) ) * ( 1.0 - zswi2 )
438               zswi1              = 1. - zswi2 * zswi12
439               zfracs             = MIN ( zswi1  * 0.12 + zswi12 * ( 0.8925 + 0.0568 * LOG( 100.0 * zgrr ) )   &
440                  &               + zswi2  * 0.26 / ( 0.26 + 0.74 * EXP ( - 724300.0 * zgrr ) )  , 0.5 )
[825]441
[4688]442               ii = MOD( npb(ji) - 1, jpi ) + 1 ; ij = ( npb(ji) - 1 ) / jpi + 1
[825]443
[4688]444               s_i_new(ji)        = zswitch_sal * zfracs * sss_m(ii,ij)  &  ! New ice salinity
[4872]445                                  + ( 1. - zswitch_sal ) * sm_i_1d(ji) 
[4688]446               ! New ice growth
[5123]447               ztmelts            = - tmut * s_i_new(ji) + rt0          ! New ice melting point (K)
[825]448
[4872]449               zt_i_new           = zswitch_sal * t_bo_1d(ji) + ( 1. - zswitch_sal) * t_i_1d(ji, nlay_i)
[4688]450               
451               zEi                = cpic * ( zt_i_new - ztmelts ) &     ! Specific enthalpy of forming ice (J/kg, <0)     
[5123]452                  &               - lfus * ( 1.0 - ( ztmelts - rt0 ) / ( zt_i_new - rt0 ) )   &
453                  &               + rcp  * ( ztmelts-rt0 )         
[4688]454
[4872]455               zEw                = rcp  * ( t_bo_1d(ji) - rt0 )         ! Specific enthalpy of seawater (J/kg, < 0)
[4688]456
457               zdE                = zEi - zEw                           ! Specific enthalpy difference (J/kg, <0)
458
459               dh_i_bott(ji)      = rdt_ice * MAX( 0._wp , zf_tt(ji) / ( zdE * rhoic ) )
460
[4872]461               q_i_1d(ji,nlay_i+1) = -zEi * rhoic                        ! New ice energy of melting (J/m3, >0)
[4688]462               
[5134]463            ENDIF
464         END DO
465      END DO
[4688]466
467      ! Contribution to Energy and Salt Fluxes
[825]468      DO ji = kideb, kiut
[4688]469         IF(  zf_tt(ji) < 0._wp  ) THEN
470            ! New ice growth
471                                   
472            zfmdt          = - rhoic * dh_i_bott(ji)             ! Mass flux x time step (kg/m2, < 0)
473
[5123]474            ztmelts        = - tmut * s_i_new(ji) + rt0          ! New ice melting point (K)
[4688]475           
[4872]476            zt_i_new       = zswitch_sal * t_bo_1d(ji) + ( 1. - zswitch_sal) * t_i_1d(ji, nlay_i)
[4688]477           
478            zEi            = cpic * ( zt_i_new - ztmelts ) &     ! Specific enthalpy of forming ice (J/kg, <0)     
[5123]479               &               - lfus * ( 1.0 - ( ztmelts - rt0 ) / ( zt_i_new - rt0 ) )   &
480               &               + rcp  * ( ztmelts-rt0 )         
[4688]481           
[4872]482            zEw            = rcp  * ( t_bo_1d(ji) - rt0 )         ! Specific enthalpy of seawater (J/kg, < 0)
[4688]483           
484            zdE            = zEi - zEw                           ! Specific enthalpy difference (J/kg, <0)
485           
486            ! Contribution to heat flux to the ocean [W.m-2], >0 
[4872]487            hfx_thd_1d(ji) = hfx_thd_1d(ji) + zfmdt * a_i_1d(ji) * zEw * r1_rdtice
[4688]488
489            ! Total heat flux used in this process [W.m-2], <0 
[4872]490            hfx_bog_1d(ji) = hfx_bog_1d(ji) - zfmdt * a_i_1d(ji) * zdE * r1_rdtice
[4688]491           
492            ! Contribution to salt flux, <0
[5167]493            sfx_bog_1d(ji) = sfx_bog_1d(ji) - rhoic * a_i_1d(ji) * dh_i_bott(ji) * s_i_new(ji) * r1_rdtice
[4688]494
495            ! Contribution to mass flux, <0
[4872]496            wfx_bog_1d(ji) =  wfx_bog_1d(ji) - rhoic * a_i_1d(ji) * dh_i_bott(ji) * r1_rdtice
[4688]497
498            ! update heat content (J.m-2) and layer thickness
[4872]499            qh_i_old(ji,nlay_i+1) = qh_i_old(ji,nlay_i+1) + dh_i_bott(ji) * q_i_1d(ji,nlay_i+1)
[4688]500            h_i_old (ji,nlay_i+1) = h_i_old (ji,nlay_i+1) + dh_i_bott(ji)
[825]501         ENDIF
502      END DO
503
[4688]504      !----------------
505      ! 4.2 Basal melt
506      !----------------
507      zdeltah(:,:) = 0._wp ! important
[825]508      DO jk = nlay_i, 1, -1
509         DO ji = kideb, kiut
[5202]510            IF(  zf_tt(ji)  >  0._wp  .AND. jk > icount(ji,jk) ) THEN   ! do not calculate where layer has already disappeared by surface melting
[825]511
[5123]512               ztmelts = - tmut * s_i_1d(ji,jk) + rt0  ! Melting point of layer jk (K)
[825]513
[4872]514               IF( t_i_1d(ji,jk) >= ztmelts ) THEN !!! Internal melting
[825]515
[5123]516                  zEi               = - q_i_1d(ji,jk) * r1_rhoic    ! Specific enthalpy of melting ice (J/kg, <0)
[4688]517                  zdE               = 0._wp                         ! Specific enthalpy difference   (J/kg, <0)
518                                                                    ! set up at 0 since no energy is needed to melt water...(it is already melted)
[5202]519                  zdeltah   (ji,jk) = MIN( 0._wp , - zh_i(ji,jk) )  ! internal melting occurs when the internal temperature is above freezing     
520                                                                    ! this should normally not happen, but sometimes, heat diffusion leads to this
[4161]521
[4688]522                  dh_i_bott (ji)    = dh_i_bott(ji) + zdeltah(ji,jk)
[825]523
[5202]524                  zfmdt             = - zdeltah(ji,jk) * rhoic      ! Mass flux x time step > 0
[825]525
[4688]526                  ! Contribution to heat flux to the ocean [W.m-2], <0 (ice enthalpy zEi is "sent" to the ocean)
[4872]527                  hfx_res_1d(ji) = hfx_res_1d(ji) + zfmdt * a_i_1d(ji) * zEi * r1_rdtice
[825]528
[4872]529                  ! Contribution to salt flux (clem: using sm_i_1d and not s_i_1d(jk) is ok)
[5167]530                  sfx_res_1d(ji) = sfx_res_1d(ji) - rhoic * a_i_1d(ji) * zdeltah(ji,jk) * sm_i_1d(ji) * r1_rdtice
[4688]531                                   
532                  ! Contribution to mass flux
[4872]533                  wfx_res_1d(ji) =  wfx_res_1d(ji) - rhoic * a_i_1d(ji) * zdeltah(ji,jk) * r1_rdtice
[825]534
[4688]535                  ! update heat content (J.m-2) and layer thickness
[4872]536                  qh_i_old(ji,jk) = qh_i_old(ji,jk) + zdeltah(ji,jk) * q_i_1d(ji,jk)
[4688]537                  h_i_old (ji,jk) = h_i_old (ji,jk) + zdeltah(ji,jk)
[825]538
[4688]539               ELSE                               !!! Basal melting
[825]540
[5202]541                  zEi             = - q_i_1d(ji,jk) * r1_rhoic ! Specific enthalpy of melting ice (J/kg, <0)
542                  zEw             = rcp * ( ztmelts - rt0 )    ! Specific enthalpy of meltwater (J/kg, <0)
543                  zdE             = zEi - zEw                  ! Specific enthalpy difference   (J/kg, <0)
[825]544
[5202]545                  zfmdt           = - zq_bo(ji) / zdE          ! Mass flux x time step (kg/m2, >0)
[825]546
[5202]547                  zdeltah(ji,jk)  = - zfmdt * r1_rhoic         ! Gross thickness change
[4688]548
[5202]549                  zdeltah(ji,jk)  = MIN( 0._wp , MAX( zdeltah(ji,jk), - zh_i(ji,jk) ) ) ! bound thickness change
[4688]550                 
[5202]551                  zq_bo(ji)       = MAX( 0._wp , zq_bo(ji) - zdeltah(ji,jk) * rhoic * zdE ) ! update available heat. MAX is necessary for roundup errors
[4688]552
[5202]553                  dh_i_bott(ji)   = dh_i_bott(ji) + zdeltah(ji,jk)    ! Update basal melt
[4688]554
[5202]555                  zfmdt           = - zdeltah(ji,jk) * rhoic          ! Mass flux x time step > 0
[4688]556
[5202]557                  zQm             = zfmdt * zEw         ! Heat exchanged with ocean
[4688]558
559                  ! Contribution to heat flux to the ocean [W.m-2], <0 
[5202]560                  hfx_thd_1d(ji)  = hfx_thd_1d(ji) + zfmdt * a_i_1d(ji) * zEw * r1_rdtice
[4688]561
[4872]562                  ! Contribution to salt flux (clem: using sm_i_1d and not s_i_1d(jk) is ok)
[5202]563                  sfx_bom_1d(ji)  = sfx_bom_1d(ji) - rhoic *  a_i_1d(ji) * zdeltah(ji,jk) * sm_i_1d(ji) * r1_rdtice
[4688]564                 
565                  ! Total heat flux used in this process [W.m-2], >0 
[5202]566                  hfx_bom_1d(ji)  = hfx_bom_1d(ji) - zfmdt * a_i_1d(ji) * zdE * r1_rdtice
[4688]567                 
568                  ! Contribution to mass flux
[5202]569                  wfx_bom_1d(ji)  =  wfx_bom_1d(ji) - rhoic * a_i_1d(ji) * zdeltah(ji,jk) * r1_rdtice
[4688]570
571                  ! update heat content (J.m-2) and layer thickness
[4872]572                  qh_i_old(ji,jk) = qh_i_old(ji,jk) + zdeltah(ji,jk) * q_i_1d(ji,jk)
[4688]573                  h_i_old (ji,jk) = h_i_old (ji,jk) + zdeltah(ji,jk)
574               ENDIF
575           
576            ENDIF
[5123]577         END DO
578      END DO
[4688]579
[825]580      !-------------------------------------------
[4688]581      ! Update temperature, energy
[825]582      !-------------------------------------------
583      DO ji = kideb, kiut
[4872]584         ht_i_1d(ji) =  MAX( 0._wp , ht_i_1d(ji) + dh_i_bott(ji) )
[4688]585      END DO 
[825]586
[4688]587      !-------------------------------------------
588      ! 5. What to do with remaining energy
589      !-------------------------------------------
590      ! If heat still available for melting and snow remains, then melt more snow
591      !-------------------------------------------
592      zdeltah(:,:) = 0._wp ! important
593      DO ji = kideb, kiut
594         zq_rema(ji)     = zq_su(ji) + zq_bo(ji) 
[5128]595         rswitch         = 1._wp - MAX( 0._wp, SIGN( 1._wp, - ht_s_1d(ji) ) )   ! =1 if snow
[5134]596         rswitch         = rswitch * MAX( 0._wp, SIGN( 1._wp, q_s_1d(ji,1) - epsi20 ) )
597         zdeltah  (ji,1) = - rswitch * zq_rema(ji) / MAX( q_s_1d(ji,1), epsi20 )
[5128]598         zdeltah  (ji,1) = MIN( 0._wp , MAX( zdeltah(ji,1) , - ht_s_1d(ji) ) ) ! bound melting
599         dh_s_tot (ji)   = dh_s_tot(ji)  + zdeltah(ji,1)
600         ht_s_1d   (ji)  = ht_s_1d(ji)   + zdeltah(ji,1)
601       
[5134]602         zq_rema(ji)     = zq_rema(ji) + zdeltah(ji,1) * q_s_1d(ji,1)                ! update available heat (J.m-2)
[5128]603         ! heat used to melt snow
[5134]604         hfx_snw_1d(ji)  = hfx_snw_1d(ji) - zdeltah(ji,1) * a_i_1d(ji) * q_s_1d(ji,1) * r1_rdtice ! W.m-2 (>0)
[5128]605         ! Contribution to mass flux
606         wfx_snw_1d(ji)  =  wfx_snw_1d(ji) - rhosn * a_i_1d(ji) * zdeltah(ji,1) * r1_rdtice
607         !   
[4688]608         ii = MOD( npb(ji) - 1, jpi ) + 1 ; ij = ( npb(ji) - 1 ) / jpi + 1
609         ! Remaining heat flux (W.m-2) is sent to the ocean heat budget
[5123]610         hfx_out(ii,ij)  = hfx_out(ii,ij) + ( zq_rema(ji) * a_i_1d(ji) ) * r1_rdtice
[825]611
[5128]612         IF( ln_icectl .AND. zq_rema(ji) < 0. .AND. lwp ) WRITE(numout,*) 'ALERTE zq_rema <0 = ', zq_rema(ji)
[4688]613      END DO
614     
[921]615      !
616      !------------------------------------------------------------------------------|
617      !  6) Snow-Ice formation                                                       |
618      !------------------------------------------------------------------------------|
[1572]619      ! When snow load excesses Archimede's limit, snow-ice interface goes down under sea-level,
620      ! flooding of seawater transforms snow into ice dh_snowice is positive for the ice
[825]621      DO ji = kideb, kiut
[1572]622         !
[4872]623         dh_snowice(ji) = MAX(  0._wp , ( rhosn * ht_s_1d(ji) + (rhoic-rau0) * ht_i_1d(ji) ) / ( rhosn+rau0-rhoic )  )
[825]624
[5167]625         ht_i_1d(ji)    = ht_i_1d(ji) + dh_snowice(ji)
626         ht_s_1d(ji)    = ht_s_1d(ji) - dh_snowice(ji)
[825]627
[4688]628         ! Salinity of snow ice
629         ii = MOD( npb(ji) - 1, jpi ) + 1 ; ij = ( npb(ji) - 1 ) / jpi + 1
[5123]630         zs_snic = zswitch_sal * sss_m(ii,ij) * ( rhoic - rhosn ) * r1_rhoic + ( 1. - zswitch_sal ) * sm_i_1d(ji)
[825]631
632         ! entrapment during snow ice formation
[5134]633         ! new salinity difference stored (to be used in limthd_sal.F90)
[5123]634         IF (  nn_icesal == 2  ) THEN
[5134]635            rswitch = MAX( 0._wp , SIGN( 1._wp , ht_i_1d(ji) - epsi20 ) )
[4161]636            ! salinity dif due to snow-ice formation
[5134]637            dsm_i_si_1d(ji) = ( zs_snic - sm_i_1d(ji) ) * dh_snowice(ji) / MAX( ht_i_1d(ji), epsi20 ) * rswitch     
[4161]638            ! salinity dif due to bottom growth
[4688]639            IF (  zf_tt(ji)  < 0._wp ) THEN
[5134]640               dsm_i_se_1d(ji) = ( s_i_new(ji) - sm_i_1d(ji) ) * dh_i_bott(ji) / MAX( ht_i_1d(ji), epsi20 ) * rswitch
[4161]641            ENDIF
642         ENDIF
[825]643
[4688]644         ! Contribution to energy flux to the ocean [J/m2], >0 (if sst<0)
645         ii = MOD( npb(ji) - 1, jpi ) + 1 ; ij = ( npb(ji) - 1 ) / jpi + 1
646         zfmdt          = ( rhosn - rhoic ) * MAX( dh_snowice(ji), 0._wp )    ! <0
647         zsstK          = sst_m(ii,ij) + rt0                               
648         zEw            = rcp * ( zsstK - rt0 )
649         zQm            = zfmdt * zEw 
650         
651         ! Contribution to heat flux
[4872]652         hfx_thd_1d(ji) = hfx_thd_1d(ji) + zfmdt * a_i_1d(ji) * zEw * r1_rdtice 
[825]653
[4688]654         ! Contribution to salt flux
[4872]655         sfx_sni_1d(ji) = sfx_sni_1d(ji) + sss_m(ii,ij) * a_i_1d(ji) * zfmdt * r1_rdtice 
[4688]656         
657         ! Contribution to mass flux
658         ! All snow is thrown in the ocean, and seawater is taken to replace the volume
[4872]659         wfx_sni_1d(ji) = wfx_sni_1d(ji) - a_i_1d(ji) * dh_snowice(ji) * rhoic * r1_rdtice
660         wfx_snw_1d(ji) = wfx_snw_1d(ji) + a_i_1d(ji) * dh_snowice(ji) * rhosn * r1_rdtice
[4688]661
662         ! update heat content (J.m-2) and layer thickness
[4872]663         qh_i_old(ji,0) = qh_i_old(ji,0) + dh_snowice(ji) * q_s_1d(ji,1) + zfmdt * zEw
[4688]664         h_i_old (ji,0) = h_i_old (ji,0) + dh_snowice(ji)
665         
[5134]666      END DO
[825]667
[2715]668      !
[4688]669      !-------------------------------------------
670      ! Update temperature, energy
671      !-------------------------------------------
672      DO ji = kideb, kiut
[4990]673         rswitch     =  1.0 - MAX( 0._wp , SIGN( 1._wp , - ht_i_1d(ji) ) ) 
[5123]674         t_su_1d(ji) =  rswitch * t_su_1d(ji) + ( 1.0 - rswitch ) * rt0
[5134]675      END DO
[4688]676
677      DO jk = 1, nlay_s
678         DO ji = kideb,kiut
679            ! mask enthalpy
[5134]680            rswitch       = 1._wp - MAX(  0._wp , SIGN( 1._wp, - ht_s_1d(ji) )  )
681            q_s_1d(ji,jk) = rswitch * q_s_1d(ji,jk)
[4872]682            ! recalculate t_s_1d from q_s_1d
[5134]683            t_s_1d(ji,jk) = rt0 + rswitch * ( - q_s_1d(ji,jk) / ( rhosn * cpic ) + lfus / cpic )
[4688]684         END DO
685      END DO
[5134]686
687      ! --- ensure that a_i = 0 where ht_i = 0 ---
688      WHERE( ht_i_1d == 0._wp ) a_i_1d = 0._wp
[5123]689     
[5167]690      CALL wrk_dealloc( jpij, zqprec, zq_su, zq_bo, zf_tt, zq_rema )
[4688]691      CALL wrk_dealloc( jpij, zdh_s_mel, zdh_s_pre, zdh_s_sub, zqh_i, zqh_s, zq_s )
[5202]692      CALL wrk_dealloc( jpij, nlay_i, zdeltah, zh_i )
693      CALL wrk_dealloc( jpij, nlay_i, icount )
[2715]694      !
[4161]695      !
[921]696   END SUBROUTINE lim_thd_dh
[1572]697   
[825]698#else
[1572]699   !!----------------------------------------------------------------------
700   !!   Default option                               NO  LIM3 sea-ice model
701   !!----------------------------------------------------------------------
[825]702CONTAINS
703   SUBROUTINE lim_thd_dh          ! Empty routine
704   END SUBROUTINE lim_thd_dh
705#endif
[1572]706
707   !!======================================================================
[921]708END MODULE limthd_dh
Note: See TracBrowser for help on using the repository browser.