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_r7881_no_wrk_alloc/NEMOGCM/NEMO/LIM_SRC_3 – NEMO

source: branches/2017/dev_r7881_no_wrk_alloc/NEMOGCM/NEMO/LIM_SRC_3/limthd_lac.F90 @ 7910

Last change on this file since 7910 was 7910, checked in by timgraham, 7 years ago

All wrk_alloc removed

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