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 branches/2016/v3_6_CMIP6_ice_diagnostics/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: branches/2016/v3_6_CMIP6_ice_diagnostics/NEMOGCM/NEMO/LIM_SRC_3/limthd_dh.F90 @ 7597

Last change on this file since 7597 was 7506, checked in by vancop, 7 years ago

Commit a first set of modifications

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