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

source: branches/2017/dev_r8183_ICEMODEL/NEMOGCM/NEMO/LIM_SRC_3/icethd_do.F90 @ 8559

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

almost useless commits

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