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_lac.F90 in branches/2017/dev_r8183_ICEMODEL/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: branches/2017/dev_r8183_ICEMODEL/NEMOGCM/NEMO/LIM_SRC_3/limthd_lac.F90 @ 8379

Last change on this file since 8379 was 8379, checked in by clem, 7 years ago

more cleaning

  • Property svn:keywords set to Id
File size: 23.4 KB
RevLine 
[825]1MODULE limthd_lac
2   !!======================================================================
3   !!                       ***  MODULE limthd_lac   ***
4   !!                lateral thermodynamic growth of the ice
5   !!======================================================================
[2715]6   !! History :  LIM  ! 2005-12 (M. Vancoppenolle)  Original code
7   !!             -   ! 2006-01 (M. Vancoppenolle)  add ITD
8   !!            3.0  ! 2007-07 (M. Vancoppenolle)  Mass and energy conservation tested
9   !!            4.0  ! 2011-02 (G. Madec) dynamical allocation
10   !!----------------------------------------------------------------------
[888]11#if defined key_lim3
[825]12   !!----------------------------------------------------------------------
[2528]13   !!   'key_lim3'                                      LIM3 sea-ice model
14   !!----------------------------------------------------------------------
[3625]15   !!   lim_lat_acr   : lateral accretion of ice
[2528]16   !!----------------------------------------------------------------------
[3625]17   USE par_oce        ! ocean parameters
18   USE dom_oce        ! domain variables
19   USE phycst         ! physical constants
[8316]20   USE sbc_oce , ONLY : sss_m
21   USE sbc_ice , ONLY : utau_ice, vtau_ice
[3625]22   USE thd_ice        ! LIM thermodynamics
23   USE ice            ! LIM variables
24   USE limtab         ! LIM 2D <==> 1D
25   USE limcons        ! LIM conservation
[8378]26   USE limthd_ent
27   USE limvar
28   !
[3625]29   USE in_out_manager ! I/O manager
30   USE lib_mpp        ! MPP library
[4833]31   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
[3625]32   USE lib_fortran    ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) 
[921]33
[825]34   IMPLICIT NONE
35   PRIVATE
36
37   PUBLIC lim_thd_lac     ! called by lim_thd
38
39   !!----------------------------------------------------------------------
[4161]40   !! NEMO/LIM3 4.0 , UCL - NEMO Consortium (2011)
[1156]41   !! $Id$
[2715]42   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
[825]43   !!----------------------------------------------------------------------
44CONTAINS
[921]45
[825]46   SUBROUTINE lim_thd_lac
47      !!-------------------------------------------------------------------
48      !!               ***   ROUTINE lim_thd_lac  ***
49      !! 
50      !! ** Purpose : Computation of the evolution of the ice thickness and
51      !!      concentration as a function of the heat balance in the leads.
52      !!      It is only used for lateral accretion
53      !!       
54      !! ** Method  : Ice is formed in the open water when ocean lose heat
55      !!      (heat budget of open water Bl is negative) .
56      !!      Computation of the increase of 1-A (ice concentration) fol-
57      !!      lowing the law :
58      !!      (dA/dt)acc = F[ (1-A)/(1-a) ] * [ Bl / (Li*h0) ]
59      !!       where - h0 is the thickness of ice created in the lead
60      !!             - a is a minimum fraction for leads
61      !!             - F is a monotonic non-increasing function defined as:
62      !!                  F(X)=( 1 - X**exld )**(1.0/exld)
63      !!             - exld is the exponent closure rate (=2 default val.)
64      !!
65      !! ** Action : - Adjustment of snow and ice thicknesses and heat
66      !!                content in brine pockets
67      !!             - Updating ice internal temperature
68      !!             - Computation of variation of ice volume and mass
[8313]69      !!             - Computation of a_i after lateral accretion and
[4872]70      !!               update ht_s_1d, ht_i_1d and tbif_1d(:,:)     
[825]71      !!------------------------------------------------------------------------
[7646]72      INTEGER  ::   ji,jj,jk,jl      ! dummy loop indices
[8327]73      INTEGER  ::   iter     !   -       -
[7646]74      REAL(wp) ::   ztmelts, zdv, zfrazb, zweight, zde                          ! local scalars
[6416]75      REAL(wp) ::   zgamafr, zvfrx, zvgx, ztaux, ztwogp, zf                     !   -      -
[2715]76      REAL(wp) ::   ztenagm, zvfry, zvgy, ztauy, zvrel2, zfp, zsqcd , zhicrit   !   -      -
[4688]77
78      REAL(wp) ::   zQm          ! enthalpy exchanged with the ocean (J/m2, >0 towards ocean)
79      REAL(wp) ::   zEi          ! sea ice specific enthalpy (J/kg)
80      REAL(wp) ::   zEw          ! seawater specific enthalpy (J/kg)
81      REAL(wp) ::   zfmdt        ! mass flux x time step (kg/m2, >0 towards ocean)
82     
83      REAL(wp) ::   zv_newfra
[8379]84      REAL(wp) ::   zvi_b, zsmv_b, zei_b, zfs_b, zfw_b, zft_b ! conservation check
[4688]85 
[8373]86      INTEGER , DIMENSION(jpij) ::   jcat        ! indexes of categories where new ice grows
87      REAL(wp), DIMENSION(jpij) ::   zswinew     ! switch for new ice or not
[825]88
[8373]89      REAL(wp), DIMENSION(jpij) ::   zv_newice   ! volume of accreted ice
90      REAL(wp), DIMENSION(jpij) ::   za_newice   ! fractional area of accreted ice
91      REAL(wp), DIMENSION(jpij) ::   zh_newice   ! thickness of accreted ice
92      REAL(wp), DIMENSION(jpij) ::   ze_newice   ! heat content of accreted ice
93      REAL(wp), DIMENSION(jpij) ::   zs_newice   ! salinity of accreted ice
94      REAL(wp), DIMENSION(jpij) ::   zo_newice   ! age of accreted ice
95      REAL(wp), DIMENSION(jpij) ::   zdv_res     ! residual volume in case of excessive heat budget
96      REAL(wp), DIMENSION(jpij) ::   zda_res     ! residual area in case of excessive heat budget
97      REAL(wp), DIMENSION(jpij) ::   zat_i_1d    ! total ice fraction   
98      REAL(wp), DIMENSION(jpij) ::   zv_frazb    ! accretion of frazil ice at the ice bottom
99      REAL(wp), DIMENSION(jpij) ::   zvrel_1d    ! relative ice / frazil velocity (1D vector)
[825]100
[8373]101      REAL(wp), DIMENSION(jpij,jpl) ::   zv_b      ! old volume of ice in category jl
102      REAL(wp), DIMENSION(jpij,jpl) ::   za_b      ! old area of ice in category jl
103      REAL(wp), DIMENSION(jpij,jpl) ::   za_i_1d   ! 1-D version of a_i
104      REAL(wp), DIMENSION(jpij,jpl) ::   zv_i_1d   ! 1-D version of v_i
105      REAL(wp), DIMENSION(jpij,jpl) ::   zsmv_i_1d ! 1-D version of smv_i
[825]106
[8373]107      REAL(wp), DIMENSION(jpij,nlay_i,jpl) ::   ze_i_1d !: 1-D version of e_i
[825]108
[8373]109      REAL(wp), DIMENSION(jpi,jpj) ::   zvrel     ! relative ice / frazil velocity
[5123]110
[6416]111      REAL(wp) :: zcai = 1.4e-3_wp                     ! ice-air drag (clem: should be dependent on coupling/forcing used)
[3294]112      !!-----------------------------------------------------------------------!
[825]113
[8379]114      IF( ln_limdiachk ) CALL lim_cons_hsm(0, 'limthd_lac', zvi_b, zsmv_b, zei_b, zfw_b, zfs_b, zft_b)
115
[5167]116      CALL lim_var_agg(1)
117      CALL lim_var_glo2eqv
[825]118
[921]119      !------------------------------------------------------------------------------!
120      ! 3) Collection thickness of ice formed in leads and polynyas
121      !------------------------------------------------------------------------------!   
[865]122      ! hicol is the thickness of new ice formed in open water
[6416]123      ! hicol can be either prescribed (frazswi = 0) or computed (frazswi = 1)
[825]124      ! Frazil ice forms in open water, is transported by wind
125      ! accumulates at the edge of the consolidated ice edge
126      ! where it forms aggregates of a specific thickness called
127      ! collection thickness.
128
[865]129      ! Note : the following algorithm currently breaks vectorization
130      !
131
[7753]132      zvrel(:,:) = 0._wp
[825]133
[7753]134      ! Default new ice thickness
135      WHERE( qlead(:,:) < 0._wp ) ; hicol(:,:) = rn_hnewice
136      ELSEWHERE                   ; hicol(:,:) = 0._wp
137      END WHERE
[825]138
[5123]139      IF( ln_frazil ) THEN
[825]140
[921]141         !--------------------
142         ! Physical constants
143         !--------------------
[7753]144         hicol(:,:) = 0._wp
[825]145
[921]146         zhicrit = 0.04 ! frazil ice thickness
147         ztwogp  = 2. * rau0 / ( grav * 0.3 * ( rau0 - rhoic ) ) ! reduced grav
[5123]148         zsqcd   = 1.0 / SQRT( 1.3 * zcai ) ! 1/SQRT(airdensity*drag)
[921]149         zgamafr = 0.03
[825]150
[7646]151         DO jj = 2, jpjm1
152            DO ji = 2, jpim1
153               IF ( qlead(ji,jj) < 0._wp .AND. tau_icebfr(ji,jj) == 0._wp ) THEN ! activated if cooling and no landfast
[921]154                  !-------------
155                  ! Wind stress
156                  !-------------
157                  ! C-grid wind stress components
[5123]158                  ztaux         = ( utau_ice(ji-1,jj  ) * umask(ji-1,jj  ,1)   &
159                     &          +   utau_ice(ji  ,jj  ) * umask(ji  ,jj  ,1) ) * 0.5_wp
160                  ztauy         = ( vtau_ice(ji  ,jj-1) * vmask(ji  ,jj-1,1)   &
161                     &          +   vtau_ice(ji  ,jj  ) * vmask(ji  ,jj  ,1) ) * 0.5_wp
[921]162                  ! Square root of wind stress
[6416]163                  ztenagm       =  SQRT( SQRT( ztaux * ztaux + ztauy * ztauy ) )
[825]164
[921]165                  !---------------------
166                  ! Frazil ice velocity
167                  !---------------------
[4990]168                  rswitch = MAX( 0._wp, SIGN( 1._wp , ztenagm - epsi10 ) )
169                  zvfrx   = rswitch * zgamafr * zsqcd * ztaux / MAX( ztenagm, epsi10 )
170                  zvfry   = rswitch * zgamafr * zsqcd * ztauy / MAX( ztenagm, epsi10 )
[825]171
[921]172                  !-------------------
173                  ! Pack ice velocity
174                  !-------------------
175                  ! C-grid ice velocity
[7646]176                  zvgx    = ( u_ice(ji-1,jj  ) * umask(ji-1,jj  ,1)  + u_ice(ji,jj) * umask(ji,jj,1) ) * 0.5_wp
177                  zvgy    = ( v_ice(ji  ,jj-1) * vmask(ji  ,jj-1,1)  + v_ice(ji,jj) * vmask(ji,jj,1) ) * 0.5_wp
[825]178
[921]179                  !-----------------------------------
180                  ! Relative frazil/pack ice velocity
181                  !-----------------------------------
182                  ! absolute relative velocity
[7646]183                  rswitch      = MAX( 0._wp, SIGN( 1._wp , at_i(ji,jj) - epsi10 ) )
184                  zvrel2       = MAX(  ( zvfrx - zvgx ) * ( zvfrx - zvgx )   &
185                     &               + ( zvfry - zvgy ) * ( zvfry - zvgy ) , 0.15 * 0.15 ) * rswitch
[6416]186                  zvrel(ji,jj) = SQRT( zvrel2 )
[825]187
[921]188                  !---------------------
189                  ! Iterative procedure
190                  !---------------------
[6416]191                  hicol(ji,jj) = zhicrit +   ( zhicrit + 0.1 )    &
192                     &                   / ( ( zhicrit + 0.1 ) * ( zhicrit + 0.1 ) -  zhicrit * zhicrit ) * ztwogp * zvrel2
[825]193
[921]194                  iter = 1
[6416]195                  DO WHILE ( iter < 20 ) 
196                     zf  = ( hicol(ji,jj) - zhicrit ) * ( hicol(ji,jj) * hicol(ji,jj) - zhicrit * zhicrit ) -   &
197                        &    hicol(ji,jj) * zhicrit * ztwogp * zvrel2
198                     zfp = ( hicol(ji,jj) - zhicrit ) * ( 3.0 * hicol(ji,jj) + zhicrit ) - zhicrit * ztwogp * zvrel2
[825]199
[7646]200                     hicol(ji,jj) = hicol(ji,jj) - zf / MAX( zfp, epsi20 )
[921]201                     iter = iter + 1
[6416]202                  END DO
[825]203
[921]204               ENDIF ! end of selection of pixels where ice forms
[825]205
[6416]206            END DO
207         END DO 
208         !
[8373]209         CALL lbc_lnk_multi( zvrel, 'T', 1., hicol, 'T', 1.  )
[825]210
211      ENDIF ! End of computation of frazil ice collection thickness
212
[921]213      !------------------------------------------------------------------------------!
214      ! 4) Identify grid points where new ice forms
215      !------------------------------------------------------------------------------!
[825]216
217      !-------------------------------------
218      ! Select points for new ice formation
219      !-------------------------------------
[7646]220      ! This occurs if open water energy budget is negative (cooling) and there is no landfast ice
[8327]221      nidx = 0 ; idxice(:) = 0
[8332]222      DO jj = 1, jpj
223         DO ji = 1, jpi
[7646]224            IF ( qlead(ji,jj)  <  0._wp .AND. tau_icebfr(ji,jj) == 0._wp ) THEN
[8327]225               nidx = nidx + 1
226               idxice( nidx ) = (jj - 1) * jpi + ji
[825]227            ENDIF
228         END DO
229      END DO
230
231      !------------------------------
232      ! Move from 2-D to 1-D vectors
233      !------------------------------
[6416]234      ! If ocean gains heat do nothing. Otherwise compute new ice formation
[825]235
[8327]236      IF ( nidx > 0 ) THEN
[825]237
[8342]238         CALL tab_2d_1d( nidx, idxice(1:nidx), zat_i_1d  (1:nidx)     , at_i           )
[921]239         DO jl = 1, jpl
[8342]240            CALL tab_2d_1d( nidx, idxice(1:nidx), za_i_1d  (1:nidx,jl), a_i  (:,:,jl)  )
241            CALL tab_2d_1d( nidx, idxice(1:nidx), zv_i_1d  (1:nidx,jl), v_i  (:,:,jl)  )
242            CALL tab_2d_1d( nidx, idxice(1:nidx), zsmv_i_1d(1:nidx,jl), smv_i(:,:,jl)  )
[921]243            DO jk = 1, nlay_i
[8342]244               CALL tab_2d_1d( nidx, idxice(1:nidx), ze_i_1d(1:nidx,jk,jl), e_i(:,:,jk,jl)   )
[5202]245            END DO
246         END DO
[825]247
[8342]248         CALL tab_2d_1d( nidx, idxice(1:nidx), qlead_1d  (1:nidx)     , qlead       )
249         CALL tab_2d_1d( nidx, idxice(1:nidx), t_bo_1d   (1:nidx)     , t_bo        )
250         CALL tab_2d_1d( nidx, idxice(1:nidx), sfx_opw_1d(1:nidx)     , sfx_opw     )
251         CALL tab_2d_1d( nidx, idxice(1:nidx), wfx_opw_1d(1:nidx)     , wfx_opw     )
252         CALL tab_2d_1d( nidx, idxice(1:nidx), hicol_1d  (1:nidx)     , hicol       )
253         CALL tab_2d_1d( nidx, idxice(1:nidx), zvrel_1d  (1:nidx)     , zvrel       )
[834]254
[8342]255         CALL tab_2d_1d( nidx, idxice(1:nidx), hfx_thd_1d(1:nidx)     , hfx_thd     )
256         CALL tab_2d_1d( nidx, idxice(1:nidx), hfx_opw_1d(1:nidx)     , hfx_opw     )
257         CALL tab_2d_1d( nidx, idxice(1:nidx), rn_amax_1d(1:nidx)     , rn_amax_2d  )
258         CALL tab_2d_1d( nidx, idxice(1:nidx), sss_1d    (1:nidx)     , sss_m       )
[4688]259
[8327]260         !------------------------------------------------------------------------------|
261         ! 2) Convert units for ice internal energy
262         !------------------------------------------------------------------------------|
263         DO jl = 1, jpl
264            DO jk = 1, nlay_i
265               DO ji = 1, nidx
266                  IF( zv_i_1d(ji,jl) > 0._wp )   ze_i_1d(ji,jk,jl) = ze_i_1d(ji,jk,jl) / zv_i_1d(ji,jl) * REAL( nlay_i )
267               END DO
268            END DO
269         END DO
[921]270         !------------------------------------------------------------------------------!
271         ! 5) Compute thickness, salinity, enthalpy, age, area and volume of new ice
272         !------------------------------------------------------------------------------!
[825]273
[4688]274         !-----------------------------------------
275         ! Keep old ice areas and volume in memory
276         !-----------------------------------------
[8327]277         zv_b(1:nidx,:) = zv_i_1d(1:nidx,:) 
278         za_b(1:nidx,:) = za_i_1d(1:nidx,:)
[6416]279
[921]280         !----------------------
281         ! Thickness of new ice
282         !----------------------
[8327]283         zh_newice(1:nidx) = hicol_1d(1:nidx)
[825]284
[921]285         !----------------------
286         ! Salinity of new ice
287         !----------------------
[5123]288         SELECT CASE ( nn_icesal )
[3625]289         CASE ( 1 )                    ! Sice = constant
[8327]290            zs_newice(1:nidx) = rn_icesal
[3625]291         CASE ( 2 )                    ! Sice = F(z,t) [Vancoppenolle et al (2005)]
[8327]292            DO ji = 1, nidx
293               zs_newice(ji) = MIN(  4.606 + 0.91 / zh_newice(ji) , rn_simax , 0.5 * sss_1d(ji) )
[3625]294            END DO
295         CASE ( 3 )                    ! Sice = F(z) [multiyear ice]
[8327]296            zs_newice(1:nidx) =   2.3
[3625]297         END SELECT
[825]298
[921]299         !-------------------------
300         ! Heat content of new ice
301         !-------------------------
302         ! We assume that new ice is formed at the seawater freezing point
[8327]303         DO ji = 1, nidx
[5123]304            ztmelts       = - tmut * zs_newice(ji) + rt0                  ! Melting point (K)
[6416]305            ze_newice(ji) =   rhoic * (  cpic * ( ztmelts - t_bo_1d(ji) )                                         &
[5123]306               &                       + lfus * ( 1.0 - ( ztmelts - rt0 ) / MIN( t_bo_1d(ji) - rt0, -epsi10 ) )   &
307               &                       - rcp  *         ( ztmelts - rt0 )  )
308         END DO
[4688]309
[921]310         !----------------
311         ! Age of new ice
312         !----------------
[8327]313         DO ji = 1, nidx
[3625]314            zo_newice(ji) = 0._wp
[5202]315         END DO
[825]316
[921]317         !-------------------
318         ! Volume of new ice
319         !-------------------
[8327]320         DO ji = 1, nidx
[825]321
[5123]322            zEi           = - ze_newice(ji) * r1_rhoic             ! specific enthalpy of forming ice [J/kg]
[4688]323
[5123]324            zEw           = rcp * ( t_bo_1d(ji) - rt0 )            ! specific enthalpy of seawater at t_bo_1d [J/kg]
[4688]325                                                                   ! clem: we suppose we are already at the freezing point (condition qlead<0 is satisfyied)
326                                                                   
327            zdE           = zEi - zEw                              ! specific enthalpy difference [J/kg]
328                                             
329            zfmdt         = - qlead_1d(ji) / zdE                   ! Fm.dt [kg/m2] (<0)
330                                                                   ! clem: we use qlead instead of zqld (limthd) because we suppose we are at the freezing point   
[5123]331            zv_newice(ji) = - zfmdt * r1_rhoic
[4688]332
333            zQm           = zfmdt * zEw                            ! heat to the ocean >0 associated with mass flux 
334
335            ! Contribution to heat flux to the ocean [W.m-2], >0 
336            hfx_thd_1d(ji) = hfx_thd_1d(ji) + zfmdt * zEw * r1_rdtice
337            ! Total heat flux used in this process [W.m-2] 
338            hfx_opw_1d(ji) = hfx_opw_1d(ji) - zfmdt * zdE * r1_rdtice
339            ! mass flux
340            wfx_opw_1d(ji) = wfx_opw_1d(ji) - zv_newice(ji) * rhoic * r1_rdtice
341            ! salt flux
342            sfx_opw_1d(ji) = sfx_opw_1d(ji) - zv_newice(ji) * rhoic * zs_newice(ji) * r1_rdtice
[6416]343         END DO
344         
345         zv_frazb(:) = 0._wp
346         IF( ln_frazil ) THEN
[921]347            ! A fraction zfrazb of frazil ice is accreted at the ice bottom
[8327]348            DO ji = 1, nidx
[6416]349               rswitch       = 1._wp - MAX( 0._wp, SIGN( 1._wp , - zat_i_1d(ji) ) )
350               zfrazb        = rswitch * ( TANH( rn_Cfrazb * ( zvrel_1d(ji) - rn_vfrazb ) ) + 1.0 ) * 0.5 * rn_maxfrazb
351               zv_frazb(ji)  =         zfrazb   * zv_newice(ji)
352               zv_newice(ji) = ( 1.0 - zfrazb ) * zv_newice(ji)
353            END DO
354         END IF
355         
[921]356         !-----------------
357         ! Area of new ice
358         !-----------------
[8327]359         DO ji = 1, nidx
[3625]360            za_newice(ji) = zv_newice(ji) / zh_newice(ji)
[4688]361         END DO
[825]362
[921]363         !------------------------------------------------------------------------------!
364         ! 6) Redistribute new ice area and volume into ice categories                  !
365         !------------------------------------------------------------------------------!
[825]366
[4688]367         !------------------------
368         ! 6.1) lateral ice growth
369         !------------------------
[921]370         ! If lateral ice growth gives an ice concentration gt 1, then
[3625]371         ! we keep the excessive volume in memory and attribute it later to bottom accretion
[8327]372         DO ji = 1, nidx
[6403]373            IF ( za_newice(ji) >  ( rn_amax_1d(ji) - zat_i_1d(ji) ) ) THEN
374               zda_res(ji)   = za_newice(ji) - ( rn_amax_1d(ji) - zat_i_1d(ji) )
[3625]375               zdv_res(ji)   = zda_res  (ji) * zh_newice(ji) 
376               za_newice(ji) = za_newice(ji) - zda_res  (ji)
377               zv_newice(ji) = zv_newice(ji) - zdv_res  (ji)
[921]378            ELSE
[3625]379               zda_res(ji) = 0._wp
380               zdv_res(ji) = 0._wp
[921]381            ENDIF
[4688]382         END DO
[825]383
[4688]384         ! find which category to fill
385         zat_i_1d(:) = 0._wp
[921]386         DO jl = 1, jpl
[8327]387            DO ji = 1, nidx
[4688]388               IF( zh_newice(ji) > hi_max(jl-1) .AND. zh_newice(ji) <= hi_max(jl) ) THEN
389                  za_i_1d (ji,jl) = za_i_1d (ji,jl) + za_newice(ji)
390                  zv_i_1d (ji,jl) = zv_i_1d (ji,jl) + zv_newice(ji)
391                  jcat    (ji)    = jl
[921]392               ENDIF
[4688]393               zat_i_1d(ji) = zat_i_1d(ji) + za_i_1d  (ji,jl)
[3625]394            END DO
395         END DO
[825]396
[4688]397         ! Heat content
[8327]398         DO ji = 1, nidx
[4688]399            jl = jcat(ji)                                                    ! categroy in which new ice is put
[4872]400            zswinew  (ji) = MAX( 0._wp , SIGN( 1._wp , - za_b(ji,jl) ) )   ! 0 if old ice
[921]401         END DO
[825]402
[921]403         DO jk = 1, nlay_i
[8327]404            DO ji = 1, nidx
[4688]405               jl = jcat(ji)
[4990]406               rswitch = MAX( 0._wp, SIGN( 1._wp , zv_i_1d(ji,jl) - epsi20 ) )
[6416]407               ze_i_1d(ji,jk,jl) = zswinew(ji)   *   ze_newice(ji) +                                                    &
[4872]408                  &        ( 1.0 - zswinew(ji) ) * ( ze_newice(ji) * zv_newice(ji) + ze_i_1d(ji,jk,jl) * zv_b(ji,jl) )  &
[4990]409                  &        * rswitch / MAX( zv_i_1d(ji,jl), epsi20 )
[2715]410            END DO
411         END DO
[825]412
[4688]413         !------------------------------------------------
414         ! 6.2) bottom ice growth + ice enthalpy remapping
415         !------------------------------------------------
416         DO jl = 1, jpl
[825]417
[4688]418            ! for remapping
[8327]419            h_i_old (1:nidx,0:nlay_i+1) = 0._wp
420            eh_i_old(1:nidx,0:nlay_i+1) = 0._wp
[921]421            DO jk = 1, nlay_i
[8327]422               DO ji = 1, nidx
[5123]423                  h_i_old (ji,jk) = zv_i_1d(ji,jl) * r1_nlay_i
[8325]424                  eh_i_old(ji,jk) = ze_i_1d(ji,jk,jl) * h_i_old(ji,jk)
[2715]425               END DO
426            END DO
[4688]427
428            ! new volumes including lateral/bottom accretion + residual
[8327]429            DO ji = 1, nidx
[4990]430               rswitch        = MAX( 0._wp, SIGN( 1._wp , zat_i_1d(ji) - epsi20 ) )
431               zv_newfra      = rswitch * ( zdv_res(ji) + zv_frazb(ji) ) * za_i_1d(ji,jl) / MAX( zat_i_1d(ji) , epsi20 )
432               za_i_1d(ji,jl) = rswitch * za_i_1d(ji,jl)               
[4688]433               zv_i_1d(ji,jl) = zv_i_1d(ji,jl) + zv_newfra
434               ! for remapping
435               h_i_old (ji,nlay_i+1) = zv_newfra
[8325]436               eh_i_old(ji,nlay_i+1) = ze_newice(ji) * zv_newfra
[4688]437            ENDDO
438            ! --- Ice enthalpy remapping --- !
[8342]439            CALL lim_thd_ent( ze_i_1d(1:nidx,:,jl) ) 
[4688]440         ENDDO
441
[921]442         !-----------------
443         ! Update salinity
444         !-----------------
[4161]445         DO jl = 1, jpl
[8327]446            DO ji = 1, nidx
[4872]447               zdv   = zv_i_1d(ji,jl) - zv_b(ji,jl)
[4688]448               zsmv_i_1d(ji,jl) = zsmv_i_1d(ji,jl) + zdv * zs_newice(ji)
449            END DO
[4161]450         END DO
451
[921]452         !------------------------------------------------------------------------------!
[8327]453         ! 8) Change units for e_i
454         !------------------------------------------------------------------------------!   
455         DO jl = 1, jpl
456            DO jk = 1, nlay_i
457               DO ji = 1, nidx
458                  ze_i_1d(ji,jk,jl) = ze_i_1d(ji,jk,jl) * zv_i_1d(ji,jl) * r1_nlay_i 
459               END DO
460            END DO
461         END DO
462         !------------------------------------------------------------------------------!
[4688]463         ! 7) Change 2D vectors to 1D vectors
[921]464         !------------------------------------------------------------------------------!
465         DO jl = 1, jpl
[8342]466            CALL tab_1d_2d( nidx, idxice(1:nidx), za_i_1d (1:nidx,jl), a_i (:,:,jl)  )
467            CALL tab_1d_2d( nidx, idxice(1:nidx), zv_i_1d (1:nidx,jl), v_i (:,:,jl)  )
468            CALL tab_1d_2d( nidx, idxice(1:nidx), zsmv_i_1d(1:nidx,jl), smv_i (:,:,jl)   )
[921]469            DO jk = 1, nlay_i
[8342]470               CALL tab_1d_2d( nidx, idxice(1:nidx), ze_i_1d(1:nidx,jk,jl), e_i(:,:,jk,jl)  )
[2715]471            END DO
472         END DO
[8342]473         CALL tab_1d_2d( nidx, idxice(1:nidx), sfx_opw_1d(1:nidx), sfx_opw  )
474         CALL tab_1d_2d( nidx, idxice(1:nidx), wfx_opw_1d(1:nidx), wfx_opw  )
475         CALL tab_1d_2d( nidx, idxice(1:nidx), hfx_thd_1d(1:nidx), hfx_thd  )
476         CALL tab_1d_2d( nidx, idxice(1:nidx), hfx_opw_1d(1:nidx), hfx_opw  )
[2715]477         !
[8327]478      ENDIF ! nidx > 0
[2715]479      !
[8379]480      IF( ln_limdiachk ) CALL lim_cons_hsm(1, 'limthd_lac', zvi_b, zsmv_b, zei_b, zfw_b, zfs_b, zft_b)
481      !
[825]482   END SUBROUTINE lim_thd_lac
483
484#else
[2715]485   !!----------------------------------------------------------------------
486   !!   Default option                               NO  LIM3 sea-ice model
487   !!----------------------------------------------------------------------
[825]488CONTAINS
489   SUBROUTINE lim_thd_lac           ! Empty routine
490   END SUBROUTINE lim_thd_lac
491#endif
[2715]492
493   !!======================================================================
[825]494END MODULE limthd_lac
Note: See TracBrowser for help on using the repository browser.