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

source: branches/UKMO/r6232_tracer_advection/NEMOGCM/NEMO/LIM_SRC_3/limthd_lac.F90 @ 9295

Last change on this file since 9295 was 9295, checked in by jcastill, 6 years ago

Remove svn keywords

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