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

source: branches/2014/dev_r4650_UKMO13_CICE_changes_take2/NEMOGCM/NEMO/LIM_SRC_3/limthd_lac.F90 @ 4921

Last change on this file since 4921 was 4921, checked in by timgraham, 9 years ago

merged with revision 4879 of trunk

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