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.F90 in NEMO/branches/2021/dev_r14383_PISCES_NEWDEV_PISCO/src/ICE – NEMO

source: NEMO/branches/2021/dev_r14383_PISCES_NEWDEV_PISCO/src/ICE/icethd.F90 @ 15349

Last change on this file since 15349 was 15349, checked in by cetlod, 3 years ago

dev_PISCO : update with trunk@15348

  • Property svn:keywords set to Id
File size: 39.8 KB
Line 
1MODULE icethd
2   !!======================================================================
3   !!                  ***  MODULE icethd   ***
4   !!   sea-ice : master routine for thermodynamics
5   !!======================================================================
6   !! History :  1.0  !  2000-01  (M.A. Morales Maqueda, H. Goosse, T. Fichefet) original code 1D
7   !!            4.0  !  2018     (many people)       SI3 [aka Sea Ice cube]
8   !!----------------------------------------------------------------------
9#if defined key_si3
10   !!----------------------------------------------------------------------
11   !!   'key_si3'                                       SI3 sea-ice model
12   !!----------------------------------------------------------------------
13   !!   ice_thd       : thermodynamics of sea ice
14   !!   ice_thd_init  : initialisation of sea-ice thermodynamics
15   !!----------------------------------------------------------------------
16   USE phycst         ! physical constants
17   USE dom_oce        ! ocean space and time domain variables
18   USE ice            ! sea-ice: variables
19!!gm list trop longue ==>>> why not passage en argument d'appel ?
20   USE sbc_oce , ONLY : sss_m, sst_m, e3t_m, utau, vtau, ssu_m, ssv_m, frq_m, sprecip, ln_cpl
21   USE sbc_ice , ONLY : qsr_oce, qns_oce, qemp_oce, qsr_ice, qns_ice, dqns_ice, evap_ice, qprec_ice, qevap_ice, &
22      &                 qml_ice, qcn_ice, qtr_ice_top, utau_ice, vtau_ice
23   USE ice1D          ! sea-ice: thermodynamics variables
24   USE icethd_zdf     ! sea-ice: vertical heat diffusion
25   USE icethd_dh      ! sea-ice: ice-snow growth and melt
26   USE icethd_da      ! sea-ice: lateral melting
27   USE icethd_sal     ! sea-ice: salinity
28   USE icethd_ent     ! sea-ice: enthalpy redistribution
29   USE icethd_do      ! sea-ice: growth in open water
30   USE icethd_pnd     ! sea-ice: melt ponds
31   USE iceitd         ! sea-ice: remapping thickness distribution
32   USE icecor         ! sea-ice: corrections
33   USE icetab         ! sea-ice: 1D <==> 2D transformation
34   USE icevar         ! sea-ice: operations
35   USE icectl         ! sea-ice: control print
36   !
37   USE in_out_manager ! I/O manager
38   USE iom            ! I/O manager library
39   USE lib_mpp        ! MPP library
40   USE lib_fortran    ! fortran utilities (glob_sum + no signed zero)
41   USE lbclnk         ! lateral boundary conditions (or mpp links)
42   USE timing         ! Timing
43
44   IMPLICIT NONE
45   PRIVATE
46
47   PUBLIC   ice_thd         ! called by limstp module
48   PUBLIC   ice_thd_init    ! called by ice_init
49
50   !!** namelist (namthd) **
51   LOGICAL ::   ln_icedH         ! activate ice thickness change from growing/melting (T) or not (F)
52   LOGICAL ::   ln_icedA         ! activate lateral melting param. (T) or not (F)
53   LOGICAL ::   ln_icedO         ! activate ice growth in open-water (T) or not (F)
54   LOGICAL ::   ln_icedS         ! activate gravity drainage and flushing (T) or not (F)
55   LOGICAL ::   ln_leadhfx       ! heat in the leads is used to melt sea-ice before warming the ocean
56
57   !! for convergence tests
58   REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) ::   ztice_cvgerr, ztice_cvgstp
59
60   !! * Substitutions
61#  include "do_loop_substitute.h90"
62   !!----------------------------------------------------------------------
63   !! NEMO/ICE 4.0 , NEMO Consortium (2018)
64   !! $Id$
65   !! Software governed by the CeCILL license (see ./LICENSE)
66   !!----------------------------------------------------------------------
67CONTAINS
68
69   SUBROUTINE ice_thd( kt )
70      !!-------------------------------------------------------------------
71      !!                ***  ROUTINE ice_thd  ***
72      !!
73      !! ** Purpose : This routine manages ice thermodynamics
74      !!
75      !! ** Action : - computation of oceanic sensible heat flux at the ice base
76      !!                              energy budget in the leads
77      !!                              net fluxes on top of ice and of ocean
78      !!             - selection of grid cells with ice
79      !!                - call ice_thd_zdf  for vertical heat diffusion
80      !!                - call ice_thd_dh   for vertical ice growth and melt
81      !!                - call ice_thd_pnd  for melt ponds
82      !!                - call ice_thd_ent  for enthalpy remapping
83      !!                - call ice_thd_sal  for ice desalination
84      !!                - call ice_thd_temp to  retrieve temperature from ice enthalpy
85      !!                - call ice_thd_mono for extra lateral ice melt if active virtual thickness distribution
86      !!                - call ice_thd_da   for lateral ice melt
87      !!             - back to the geographic grid
88      !!                - call ice_thd_rem  for remapping thickness distribution
89      !!                - call ice_thd_do   for ice growth in leads
90      !!-------------------------------------------------------------------
91      INTEGER, INTENT(in) :: kt    ! number of iteration
92      !
93      INTEGER  :: ji, jj, jk, jl   ! dummy loop indices
94      REAL(wp) :: zfric_u, zqld, zqfr, zqfr_neg, zqfr_pos
95      REAL(wp), PARAMETER :: zfric_umin = 0._wp       ! lower bound for the friction velocity (cice value=5.e-04)
96      REAL(wp), PARAMETER :: zch        = 0.0057_wp   ! heat transfer coefficient
97      REAL(wp), DIMENSION(jpi,jpj) ::   zu_io, zv_io, zfric, zvel   ! ice-ocean velocity (m/s) and frictional velocity (m2/s2)
98      !
99      ! for collection thickness
100      INTEGER  ::   iter             !   -       -
101      REAL(wp) ::   zvfrx, zvgx, ztaux, zf                     !   -      -
102      REAL(wp) ::   ztenagm, zvfry, zvgy, ztauy, zvrel2, zfp, ztwogp   !   -      -
103      REAL(wp), PARAMETER ::   zcai = 1.4e-3_wp               ! ice-air drag (clem: should be dependent on coupling/forcing used)
104      REAL(wp), PARAMETER ::   zhicrit = 0.04_wp                                             ! frazil ice thickness
105      REAL(wp), PARAMETER ::   zsqcd   = 1.0_wp / SQRT( 1.3_wp * zcai )                      ! 1/SQRT(airdensity*drag)
106      REAL(wp), PARAMETER ::   zgamafr = 0.03_wp
107      !!-------------------------------------------------------------------
108      ! controls
109      IF( ln_timing    )   CALL timing_start('icethd')                                                             ! timing
110      IF( ln_icediachk )   CALL ice_cons_hsm(0, 'icethd', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft) ! conservation
111      IF( ln_icediachk )   CALL ice_cons2D  (0, 'icethd',  diag_v,  diag_s,  diag_t,  diag_fv,  diag_fs,  diag_ft) ! conservation
112
113      IF( kt == nit000 .AND. lwp ) THEN
114         WRITE(numout,*)
115         WRITE(numout,*) 'ice_thd: sea-ice thermodynamics'
116         WRITE(numout,*) '~~~~~~~'
117      ENDIF
118
119      ! convergence tests
120      IF( ln_zdf_chkcvg ) THEN
121         ALLOCATE( ztice_cvgerr(jpi,jpj,jpl) , ztice_cvgstp(jpi,jpj,jpl) )
122         ztice_cvgerr = 0._wp ; ztice_cvgstp = 0._wp
123      ENDIF
124
125      !---------------------------------------------!
126      ! computation of friction velocity at T points
127      !---------------------------------------------!
128      IF( ln_icedyn ) THEN
129         zu_io(:,:) = u_ice(:,:) - ssu_m(:,:)
130         zv_io(:,:) = v_ice(:,:) - ssv_m(:,:)
131         DO_2D( 0, 0, 0, 0 )
132            zfric(ji,jj) = rn_cio * ( 0.5_wp *  &
133               &                    (  zu_io(ji,jj) * zu_io(ji,jj) + zu_io(ji-1,jj) * zu_io(ji-1,jj)   &
134               &                     + zv_io(ji,jj) * zv_io(ji,jj) + zv_io(ji,jj-1) * zv_io(ji,jj-1) ) ) * tmask(ji,jj,1)
135            zvel(ji,jj) = 0.5_wp * SQRT( ( u_ice(ji-1,jj) + u_ice(ji,jj) ) * ( u_ice(ji-1,jj) + u_ice(ji,jj) ) + &
136               &                         ( v_ice(ji,jj-1) + v_ice(ji,jj) ) * ( v_ice(ji,jj-1) + v_ice(ji,jj) ) )
137         END_2D
138      ELSE      !  if no ice dynamics => transfer directly the atmospheric stress to the ocean
139         DO_2D( 0, 0, 0, 0 )
140            zfric(ji,jj) = r1_rho0 * SQRT( 0.5_wp *  &
141               &                         (  utau(ji,jj) * utau(ji,jj) + utau(ji-1,jj) * utau(ji-1,jj)   &
142               &                          + vtau(ji,jj) * vtau(ji,jj) + vtau(ji,jj-1) * vtau(ji,jj-1) ) ) * tmask(ji,jj,1)
143            zvel(ji,jj) = 0._wp
144         END_2D
145      ENDIF
146      CALL lbc_lnk( 'icethd', zfric, 'T',  1.0_wp, zvel, 'T', 1.0_wp )
147      !
148      !--------------------------------------------------------------------!
149      ! Partial computation of forcing for the thermodynamic sea ice model
150      !--------------------------------------------------------------------!
151      DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )   ! needed for qlead
152         rswitch  = tmask(ji,jj,1) * MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - epsi10 ) ) ! 0 if no ice
153         !
154         ! --- Energy received in the lead from atm-oce exchanges, zqld is defined everywhere (J.m-2) --- !
155         zqld =  tmask(ji,jj,1) * rDt_ice *  &
156            &    ( ( 1._wp - at_i_b(ji,jj) ) * qsr_oce(ji,jj) * frq_m(ji,jj) +  &
157            &      ( 1._wp - at_i_b(ji,jj) ) * qns_oce(ji,jj) + qemp_oce(ji,jj) )
158
159         ! --- Energy needed to bring ocean surface layer until its freezing, zqfr is defined everywhere (J.m-2) --- !
160         !     (mostly<0 but >0 if supercooling)
161         zqfr     = rho0 * rcp * e3t_m(ji,jj) * ( t_bo(ji,jj) - ( sst_m(ji,jj) + rt0 ) ) * tmask(ji,jj,1)  ! both < 0 (t_bo < sst) and > 0 (t_bo > sst)
162         zqfr_neg = MIN( zqfr , 0._wp )                                                                    ! only < 0
163         zqfr_pos = MAX( zqfr , 0._wp )                                                                    ! only > 0
164
165         ! --- Sensible ocean-to-ice heat flux (W/m2) --- !
166         !     (mostly>0 but <0 if supercooling)
167         zfric_u            = MAX( SQRT( zfric(ji,jj) ), zfric_umin )
168         qsb_ice_bot(ji,jj) = rswitch * rho0 * rcp * zch * zfric_u * ( ( sst_m(ji,jj) + rt0 ) - t_bo(ji,jj) )
169
170         ! upper bound for qsb_ice_bot: the heat retrieved from the ocean must be smaller than the heat necessary to reach
171         !                              the freezing point, so that we do not have SST < T_freeze
172         !                              This implies: qsb_ice_bot(ji,jj) * at_i(ji,jj) * rtdice <= - zqfr_neg
173         !                              The following formulation is ok for both normal conditions and supercooling
174         qsb_ice_bot(ji,jj) = rswitch * MIN( qsb_ice_bot(ji,jj), - zqfr_neg * r1_Dt_ice / MAX( at_i(ji,jj), epsi10 ) )
175
176         ! If conditions are always supercooled (such as at the mouth of ice-shelves), then ice grows continuously
177         ! ==> stop ice formation by artificially setting up the turbulent fluxes to 0 when volume > 20m (arbitrary)
178         IF( ( t_bo(ji,jj) - ( sst_m(ji,jj) + rt0 ) ) > 0._wp .AND. vt_i(ji,jj) >= 20._wp ) THEN
179            zqfr               = 0._wp
180            zqfr_pos           = 0._wp
181            qsb_ice_bot(ji,jj) = 0._wp
182         ENDIF
183         !
184         ! --- Energy Budget of the leads (qlead, J.m-2) --- !
185         !     qlead is the energy received from the atm. in the leads.
186         !     If warming (zqld >= 0), then the energy in the leads is used to melt ice (bottom melting) => fhld  (W/m2)
187         !     If cooling (zqld <  0), then the energy in the leads is used to grow ice in open water    => qlead (J.m-2)
188         IF( zqld >= 0._wp .AND. at_i(ji,jj) > 0._wp ) THEN
189            ! upper bound for fhld: fhld should be equal to zqld
190            !                        but we have to make sure that this heat will not make the sst drop below the freezing point
191            !                        so the max heat that can be pulled out of the ocean is zqld - qsb - zqfr_pos
192            !                        The following formulation is ok for both normal conditions and supercooling
193            fhld (ji,jj) = rswitch * MAX( 0._wp, ( zqld - zqfr_pos ) * r1_Dt_ice / MAX( at_i(ji,jj), epsi10 ) &  ! divided by at_i since this is (re)multiplied by a_i in icethd_dh.F90
194               &                                 - qsb_ice_bot(ji,jj) )
195            qlead(ji,jj) = 0._wp
196         ELSE
197            fhld (ji,jj) = 0._wp
198            ! upper bound for qlead: qlead should be equal to zqld
199            !                        but before using this heat for ice formation, we suppose that the ocean cools down till the freezing point.
200            !                        The energy for this cooling down is zqfr. Also some heat will be removed from the ocean from turbulent fluxes (qsb)
201            !                        and freezing point is reached if zqfr = zqld - qsb*a/dt
202            !                        so the max heat that can be pulled out of the ocean is zqld - qsb - zqfr
203            !                        The following formulation is ok for both normal conditions and supercooling
204            qlead(ji,jj) = MIN( 0._wp , zqld - ( qsb_ice_bot(ji,jj) * at_i(ji,jj) * rDt_ice ) - zqfr )
205         ENDIF
206         !
207         ! If ice is landfast and ice concentration reaches its max
208         ! => stop ice formation in open water
209         IF(  zvel(ji,jj) <= 5.e-04_wp .AND. at_i(ji,jj) >= rn_amax_2d(ji,jj)-epsi06 )   qlead(ji,jj) = 0._wp
210         !
211         ! If the grid cell is almost fully covered by ice (no leads)
212         ! => stop ice formation in open water
213         IF( at_i(ji,jj) >= (1._wp - epsi10) )   qlead(ji,jj) = 0._wp
214         !
215         ! If ln_leadhfx is false
216         ! => do not use energy of the leads to melt sea-ice
217         IF( .NOT.ln_leadhfx )   fhld(ji,jj) = 0._wp
218         !
219      END_2D
220
221      ! In case we bypass open-water ice formation
222      IF( .NOT. ln_icedO )  qlead(:,:) = 0._wp
223      ! In case we bypass growing/melting from top and bottom
224      IF( .NOT. ln_icedH ) THEN
225         qsb_ice_bot(:,:) = 0._wp
226         fhld       (:,:) = 0._wp
227      ENDIF
228
229
230      !---------------------------------------------------------!
231      ! Collection thickness of ice formed in leads and polynyas
232      !---------------------------------------------------------!   
233      ! ht_i_new is the thickness of new ice formed in open water
234      ! ht_i_new can be either prescribed (ln_frazil=F) or computed (ln_frazil=T)
235      ! Frazil ice forms in open water, is transported by wind, accumulates at the edge of the consolidated ice edge
236      ! where it forms aggregates of a specific thickness called collection thickness.
237      !
238      fraz_frac(:,:) = 0._wp
239      !
240      ! Default new ice thickness
241      WHERE( qlead(:,:) < 0._wp ) ! cooling
242         ht_i_new(:,:) = rn_hinew
243      ELSEWHERE
244         ht_i_new(:,:) = 0._wp
245      END WHERE
246
247      IF( ln_frazil ) THEN
248         ztwogp  = 2._wp * rho0 / ( grav * 0.3_wp * ( rho0 - rhoi ) )  ! reduced grav
249         !
250         DO_2D( 0, 0, 0, 0 )
251            IF ( qlead(ji,jj) < 0._wp ) THEN ! cooling
252               ! -- Wind stress -- !
253               ztaux         = ( utau_ice(ji-1,jj  ) * umask(ji-1,jj  ,1)   &
254                  &          +   utau_ice(ji  ,jj  ) * umask(ji  ,jj  ,1) ) * 0.5_wp
255               ztauy         = ( vtau_ice(ji  ,jj-1) * vmask(ji  ,jj-1,1)   &
256                  &          +   vtau_ice(ji  ,jj  ) * vmask(ji  ,jj  ,1) ) * 0.5_wp
257               ! Square root of wind stress
258               ztenagm       =  SQRT( SQRT( ztaux * ztaux + ztauy * ztauy ) )
259
260               ! -- Frazil ice velocity -- !
261               rswitch = MAX( 0._wp, SIGN( 1._wp , ztenagm - epsi10 ) )
262               zvfrx   = rswitch * zgamafr * zsqcd * ztaux / MAX( ztenagm, epsi10 )
263               zvfry   = rswitch * zgamafr * zsqcd * ztauy / MAX( ztenagm, epsi10 )
264
265               ! -- Pack ice velocity -- !
266               zvgx    = ( u_ice(ji-1,jj  ) * umask(ji-1,jj  ,1)  + u_ice(ji,jj) * umask(ji,jj,1) ) * 0.5_wp
267               zvgy    = ( v_ice(ji  ,jj-1) * vmask(ji  ,jj-1,1)  + v_ice(ji,jj) * vmask(ji,jj,1) ) * 0.5_wp
268
269               ! -- Relative frazil/pack ice velocity -- !
270               rswitch      = MAX( 0._wp, SIGN( 1._wp , at_i(ji,jj) - epsi10 ) )
271               zvrel2       = MAX(  ( zvfrx - zvgx ) * ( zvfrx - zvgx )   &
272                  &               + ( zvfry - zvgy ) * ( zvfry - zvgy ) , 0.15_wp * 0.15_wp ) * rswitch
273
274               !!clem
275               fraz_frac(ji,jj) = rswitch * ( TANH( rn_Cfraz * ( SQRT(zvrel2) - rn_vfraz ) ) + 1._wp ) * 0.5_wp * rn_maxfraz
276               !!clem
277               
278               ! -- new ice thickness (iterative loop) -- !
279               ht_i_new(ji,jj) = zhicrit +   ( zhicrit + 0.1_wp )    &
280                  &                      / ( ( zhicrit + 0.1_wp ) * ( zhicrit + 0.1_wp ) -  zhicrit * zhicrit ) * ztwogp * zvrel2
281               iter = 1
282               DO WHILE ( iter < 20 ) 
283                  zf  = ( ht_i_new(ji,jj) - zhicrit ) * ( ht_i_new(ji,jj) * ht_i_new(ji,jj) - zhicrit * zhicrit ) -   &
284                     &    ht_i_new(ji,jj) * zhicrit * ztwogp * zvrel2
285                  zfp = ( ht_i_new(ji,jj) - zhicrit ) * ( 3.0_wp * ht_i_new(ji,jj) + zhicrit ) - zhicrit * ztwogp * zvrel2
286
287                  ht_i_new(ji,jj) = ht_i_new(ji,jj) - zf / MAX( zfp, epsi20 )
288                  iter = iter + 1
289               END DO
290               !
291               ! bound ht_i_new (though I don't see why it should be necessary)
292               ht_i_new(ji,jj) = MAX( 0.01_wp, MIN( ht_i_new(ji,jj), hi_max(jpl) ) )
293               !
294            ELSE
295               ht_i_new(ji,jj) = 0._wp
296            ENDIF
297            !
298         END_2D
299         !
300         CALL lbc_lnk( 'icethd', fraz_frac, 'T', 1.0_wp, ht_i_new, 'T', 1.0_wp  )
301
302      ENDIF
303
304      !-------------------------------------------------------------------------------------------!
305      ! Thermodynamic computation (only on grid points covered by ice) => loop over ice categories
306      !-------------------------------------------------------------------------------------------!
307      DO jl = 1, jpl
308
309         ! select ice covered grid points
310         npti = 0 ; nptidx(:) = 0
311         DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )
312            IF ( a_i(ji,jj,jl) > epsi10 ) THEN
313               npti         = npti  + 1
314               nptidx(npti) = (jj - 1) * jpi + ji
315            ENDIF
316         END_2D
317
318         IF( npti > 0 ) THEN  ! If there is no ice, do nothing.
319            !
320                              CALL ice_thd_1d2d( jl, 1 )            ! --- Move to 1D arrays --- !
321            !                                                       ! --- & Change units of e_i, e_s from J/m2 to J/m3 --- !
322            !
323            s_i_new   (1:npti) = 0._wp ; dh_s_tot(1:npti) = 0._wp   ! --- some init --- !  (important to have them here)
324            dh_i_sum  (1:npti) = 0._wp ; dh_i_bom(1:npti) = 0._wp ; dh_i_itm  (1:npti) = 0._wp
325            dh_i_sub  (1:npti) = 0._wp ; dh_i_bog(1:npti) = 0._wp
326            dh_snowice(1:npti) = 0._wp ; dh_s_mlt(1:npti) = 0._wp
327            !
328                              CALL ice_thd_zdf                      ! --- Ice-Snow temperature --- !
329            !
330            IF( ln_icedH ) THEN                                     ! --- Growing/Melting --- !
331                              CALL ice_thd_dh                           ! Ice-Snow thickness
332                              CALL ice_thd_ent( e_i_1d(1:npti,:) )      ! Ice enthalpy remapping
333            ENDIF
334                              CALL ice_thd_sal( ln_icedS )          ! --- Ice salinity --- !
335            !
336                              CALL ice_thd_temp                     ! --- Temperature update --- !
337            !
338            IF( ln_icedH .AND. ln_virtual_itd ) &
339               &              CALL ice_thd_mono                     ! --- Extra lateral melting if virtual_itd --- !
340            !
341            IF( ln_icedA )    CALL ice_thd_da                       ! --- Lateral melting --- !
342            !
343                              CALL ice_thd_1d2d( jl, 2 )            ! --- Change units of e_i, e_s from J/m3 to J/m2 --- !
344            !                                                       ! --- & Move to 2D arrays --- !
345         ENDIF
346         !
347      END DO
348      !
349      IF( ln_icediachk )   CALL ice_cons_hsm(1, 'icethd', rdiag_v, rdiag_s, rdiag_t, rdiag_fv, rdiag_fs, rdiag_ft)
350      IF( ln_icediachk )   CALL ice_cons2D  (1, 'icethd',  diag_v,  diag_s,  diag_t,  diag_fv,  diag_fs,  diag_ft)
351      !
352      IF ( ln_pnd .AND. ln_icedH ) &
353         &                    CALL ice_thd_pnd                      ! --- Melt ponds --- !
354      !
355      IF( jpl > 1  )          CALL ice_itd_rem( kt )                ! --- Transport ice between thickness categories --- !
356      !
357      IF( ln_icedO )          CALL ice_thd_do                       ! --- Frazil ice growth in leads --- !
358      !
359                              CALL ice_cor( kt , 2 )                ! --- Corrections --- !
360      !
361      oa_i(:,:,:) = oa_i(:,:,:) + a_i(:,:,:) * rDt_ice              ! --- Ice natural aging incrementation
362      !
363      DO_2D( 0, 0, 0, 0 )                                           ! --- Ice velocity corrections
364         IF( at_i(ji,jj) == 0._wp ) THEN   ! if ice has melted
365            IF( at_i(ji+1,jj) == 0._wp )   u_ice(ji  ,jj) = 0._wp   ! right side
366            IF( at_i(ji-1,jj) == 0._wp )   u_ice(ji-1,jj) = 0._wp   ! left side
367            IF( at_i(ji,jj+1) == 0._wp )   v_ice(ji,jj  ) = 0._wp   ! upper side
368            IF( at_i(ji,jj-1) == 0._wp )   v_ice(ji,jj-1) = 0._wp   ! bottom side
369         ENDIF
370      END_2D
371      CALL lbc_lnk( 'icecor', u_ice, 'U', -1.0_wp, v_ice, 'V', -1.0_wp )
372      !
373      ! convergence tests
374      IF( ln_zdf_chkcvg ) THEN
375         CALL iom_put( 'tice_cvgerr', ztice_cvgerr ) ; DEALLOCATE( ztice_cvgerr )
376         CALL iom_put( 'tice_cvgstp', ztice_cvgstp ) ; DEALLOCATE( ztice_cvgstp )
377      ENDIF
378      !
379      ! controls
380      IF( ln_icectl )   CALL ice_prt    (kt, iiceprt, jiceprt, 1, ' - ice thermodyn. - ') ! prints
381      IF( sn_cfctl%l_prtctl )   &
382        &               CALL ice_prt3D  ('icethd')                                        ! prints
383      IF( ln_timing )   CALL timing_stop('icethd')                                        ! timing
384      !
385   END SUBROUTINE ice_thd
386
387
388   SUBROUTINE ice_thd_temp
389      !!-----------------------------------------------------------------------
390      !!                   ***  ROUTINE ice_thd_temp ***
391      !!
392      !! ** Purpose :   Computes sea ice temperature (Kelvin) from enthalpy
393      !!
394      !! ** Method  :   Formula (Bitz and Lipscomb, 1999)
395      !!-------------------------------------------------------------------
396      INTEGER  ::   ji, jk   ! dummy loop indices
397      REAL(wp) ::   ztmelts, zbbb, zccc  ! local scalar
398      !!-------------------------------------------------------------------
399      ! Recover ice temperature
400      DO jk = 1, nlay_i
401         DO ji = 1, npti
402            ztmelts       = -rTmlt * sz_i_1d(ji,jk)
403            ! Conversion q(S,T) -> T (second order equation)
404            zbbb          = ( rcp - rcpi ) * ztmelts + e_i_1d(ji,jk) * r1_rhoi - rLfus
405            zccc          = SQRT( MAX( zbbb * zbbb - 4._wp * rcpi * rLfus * ztmelts, 0._wp ) )
406            t_i_1d(ji,jk) = rt0 - ( zbbb + zccc ) * 0.5_wp * r1_rcpi
407
408            ! mask temperature
409            rswitch       = 1._wp - MAX( 0._wp , SIGN( 1._wp , - h_i_1d(ji) ) )
410            t_i_1d(ji,jk) = rswitch * t_i_1d(ji,jk) + ( 1._wp - rswitch ) * rt0
411         END DO
412      END DO
413      !
414   END SUBROUTINE ice_thd_temp
415
416
417   SUBROUTINE ice_thd_mono
418      !!-----------------------------------------------------------------------
419      !!                   ***  ROUTINE ice_thd_mono ***
420      !!
421      !! ** Purpose :   Lateral melting in case virtual_itd
422      !!                          ( dA = A/2h dh )
423      !!-----------------------------------------------------------------------
424      INTEGER  ::   ji                 ! dummy loop indices
425      REAL(wp) ::   zhi_bef            ! ice thickness before thermo
426      REAL(wp) ::   zdh_mel, zda_mel   ! net melting
427      REAL(wp) ::   zvi, zvs           ! ice/snow volumes
428      !!-----------------------------------------------------------------------
429      !
430      DO ji = 1, npti
431         zdh_mel = MIN( 0._wp, dh_i_itm(ji) + dh_i_sum(ji) + dh_i_bom(ji) + dh_snowice(ji) + dh_i_sub(ji) )
432         IF( zdh_mel < 0._wp .AND. a_i_1d(ji) > 0._wp )  THEN
433            zvi          = a_i_1d(ji) * h_i_1d(ji)
434            zvs          = a_i_1d(ji) * h_s_1d(ji)
435            ! lateral melting = concentration change
436            zhi_bef     = h_i_1d(ji) - zdh_mel
437            rswitch     = MAX( 0._wp , SIGN( 1._wp , zhi_bef - epsi20 ) )
438            zda_mel     = rswitch * a_i_1d(ji) * zdh_mel / ( 2._wp * MAX( zhi_bef, epsi20 ) )
439            a_i_1d(ji)  = MAX( epsi20, a_i_1d(ji) + zda_mel )
440            ! adjust thickness
441            h_i_1d(ji) = zvi / a_i_1d(ji)
442            h_s_1d(ji) = zvs / a_i_1d(ji)
443            ! retrieve total concentration
444            at_i_1d(ji) = a_i_1d(ji)
445         END IF
446      END DO
447      !
448   END SUBROUTINE ice_thd_mono
449
450
451   SUBROUTINE ice_thd_1d2d( kl, kn )
452      !!-----------------------------------------------------------------------
453      !!                   ***  ROUTINE ice_thd_1d2d ***
454      !!
455      !! ** Purpose :   move arrays from 1d to 2d and the reverse
456      !!-----------------------------------------------------------------------
457      INTEGER, INTENT(in) ::   kl   ! index of the ice category
458      INTEGER, INTENT(in) ::   kn   ! 1= from 2D to 1D   ;   2= from 1D to 2D
459      !
460      INTEGER ::   jk   ! dummy loop indices
461      !!-----------------------------------------------------------------------
462      !
463      SELECT CASE( kn )
464      !                    !---------------------!
465      CASE( 1 )            !==  from 2D to 1D  ==!
466         !                 !---------------------!
467         CALL tab_2d_1d( npti, nptidx(1:npti), at_i_1d(1:npti), at_i             )
468         CALL tab_2d_1d( npti, nptidx(1:npti), a_i_1d (1:npti), a_i (:,:,kl)     )
469         CALL tab_2d_1d( npti, nptidx(1:npti), h_i_1d (1:npti), h_i (:,:,kl)     )
470         CALL tab_2d_1d( npti, nptidx(1:npti), h_s_1d (1:npti), h_s (:,:,kl)     )
471         CALL tab_2d_1d( npti, nptidx(1:npti), t_su_1d(1:npti), t_su(:,:,kl)     )
472         CALL tab_2d_1d( npti, nptidx(1:npti), s_i_1d (1:npti), s_i (:,:,kl)     )
473         DO jk = 1, nlay_s
474            CALL tab_2d_1d( npti, nptidx(1:npti), t_s_1d(1:npti,jk), t_s(:,:,jk,kl)    )
475            CALL tab_2d_1d( npti, nptidx(1:npti), e_s_1d(1:npti,jk), e_s(:,:,jk,kl)    )
476         END DO
477         DO jk = 1, nlay_i
478            CALL tab_2d_1d( npti, nptidx(1:npti), t_i_1d (1:npti,jk), t_i (:,:,jk,kl)  )
479            CALL tab_2d_1d( npti, nptidx(1:npti), e_i_1d (1:npti,jk), e_i (:,:,jk,kl)  )
480            CALL tab_2d_1d( npti, nptidx(1:npti), sz_i_1d(1:npti,jk), sz_i(:,:,jk,kl)  )
481         END DO
482         !
483         CALL tab_2d_1d( npti, nptidx(1:npti), qprec_ice_1d  (1:npti), qprec_ice            )
484         CALL tab_2d_1d( npti, nptidx(1:npti), qsr_ice_1d    (1:npti), qsr_ice (:,:,kl)     )
485         CALL tab_2d_1d( npti, nptidx(1:npti), qns_ice_1d    (1:npti), qns_ice (:,:,kl)     )
486         CALL tab_2d_1d( npti, nptidx(1:npti), evap_ice_1d   (1:npti), evap_ice(:,:,kl)     )
487         CALL tab_2d_1d( npti, nptidx(1:npti), dqns_ice_1d   (1:npti), dqns_ice(:,:,kl)     )
488         CALL tab_2d_1d( npti, nptidx(1:npti), t_bo_1d       (1:npti), t_bo                 )
489         CALL tab_2d_1d( npti, nptidx(1:npti), sprecip_1d    (1:npti), sprecip              )
490         CALL tab_2d_1d( npti, nptidx(1:npti), qsb_ice_bot_1d(1:npti), qsb_ice_bot          )
491         CALL tab_2d_1d( npti, nptidx(1:npti), fhld_1d       (1:npti), fhld                 )
492
493         CALL tab_2d_1d( npti, nptidx(1:npti), qml_ice_1d    (1:npti), qml_ice    (:,:,kl) )
494         CALL tab_2d_1d( npti, nptidx(1:npti), qcn_ice_1d    (1:npti), qcn_ice    (:,:,kl) )
495         CALL tab_2d_1d( npti, nptidx(1:npti), qtr_ice_top_1d(1:npti), qtr_ice_top(:,:,kl) )
496         !
497         CALL tab_2d_1d( npti, nptidx(1:npti), wfx_snw_sni_1d(1:npti), wfx_snw_sni   )
498         CALL tab_2d_1d( npti, nptidx(1:npti), wfx_snw_sum_1d(1:npti), wfx_snw_sum   )
499         CALL tab_2d_1d( npti, nptidx(1:npti), wfx_sub_1d    (1:npti), wfx_sub       )
500         CALL tab_2d_1d( npti, nptidx(1:npti), wfx_snw_sub_1d(1:npti), wfx_snw_sub   )
501         CALL tab_2d_1d( npti, nptidx(1:npti), wfx_ice_sub_1d(1:npti), wfx_ice_sub   )
502         CALL tab_2d_1d( npti, nptidx(1:npti), wfx_err_sub_1d(1:npti), wfx_err_sub   )
503         !
504         CALL tab_2d_1d( npti, nptidx(1:npti), wfx_bog_1d (1:npti), wfx_bog          )
505         CALL tab_2d_1d( npti, nptidx(1:npti), wfx_bom_1d (1:npti), wfx_bom          )
506         CALL tab_2d_1d( npti, nptidx(1:npti), wfx_sum_1d (1:npti), wfx_sum          )
507         CALL tab_2d_1d( npti, nptidx(1:npti), wfx_sni_1d (1:npti), wfx_sni          )
508         CALL tab_2d_1d( npti, nptidx(1:npti), wfx_res_1d (1:npti), wfx_res          )
509         CALL tab_2d_1d( npti, nptidx(1:npti), wfx_spr_1d (1:npti), wfx_spr          )
510         CALL tab_2d_1d( npti, nptidx(1:npti), wfx_lam_1d (1:npti), wfx_lam          )
511         !
512         CALL tab_2d_1d( npti, nptidx(1:npti), sfx_bog_1d (1:npti), sfx_bog          )
513         CALL tab_2d_1d( npti, nptidx(1:npti), sfx_bom_1d (1:npti), sfx_bom          )
514         CALL tab_2d_1d( npti, nptidx(1:npti), sfx_sum_1d (1:npti), sfx_sum          )
515         CALL tab_2d_1d( npti, nptidx(1:npti), sfx_sni_1d (1:npti), sfx_sni          )
516         CALL tab_2d_1d( npti, nptidx(1:npti), sfx_bri_1d (1:npti), sfx_bri          )
517         CALL tab_2d_1d( npti, nptidx(1:npti), sfx_res_1d (1:npti), sfx_res          )
518         CALL tab_2d_1d( npti, nptidx(1:npti), sfx_sub_1d (1:npti), sfx_sub          )
519         CALL tab_2d_1d( npti, nptidx(1:npti), sfx_lam_1d (1:npti), sfx_lam          )
520         !
521         CALL tab_2d_1d( npti, nptidx(1:npti), hfx_thd_1d    (1:npti), hfx_thd       )
522         CALL tab_2d_1d( npti, nptidx(1:npti), hfx_spr_1d    (1:npti), hfx_spr       )
523         CALL tab_2d_1d( npti, nptidx(1:npti), hfx_sum_1d    (1:npti), hfx_sum       )
524         CALL tab_2d_1d( npti, nptidx(1:npti), hfx_bom_1d    (1:npti), hfx_bom       )
525         CALL tab_2d_1d( npti, nptidx(1:npti), hfx_bog_1d    (1:npti), hfx_bog       )
526         CALL tab_2d_1d( npti, nptidx(1:npti), hfx_dif_1d    (1:npti), hfx_dif       )
527         CALL tab_2d_1d( npti, nptidx(1:npti), hfx_opw_1d    (1:npti), hfx_opw       )
528         CALL tab_2d_1d( npti, nptidx(1:npti), hfx_snw_1d    (1:npti), hfx_snw       )
529         CALL tab_2d_1d( npti, nptidx(1:npti), hfx_sub_1d    (1:npti), hfx_sub       )
530         CALL tab_2d_1d( npti, nptidx(1:npti), hfx_res_1d    (1:npti), hfx_res       )
531         CALL tab_2d_1d( npti, nptidx(1:npti), hfx_err_dif_1d(1:npti), hfx_err_dif   )
532         !
533         ! ocean surface fields
534         CALL tab_2d_1d( npti, nptidx(1:npti), sst_1d(1:npti), sst_m )
535         CALL tab_2d_1d( npti, nptidx(1:npti), sss_1d(1:npti), sss_m )
536         CALL tab_2d_1d( npti, nptidx(1:npti), frq_m_1d(1:npti), frq_m )
537         !
538         ! to update ice age
539         CALL tab_2d_1d( npti, nptidx(1:npti), o_i_1d (1:npti), o_i (:,:,kl) )
540         CALL tab_2d_1d( npti, nptidx(1:npti), oa_i_1d(1:npti), oa_i(:,:,kl) )
541         !
542         ! --- Change units of e_i, e_s from J/m2 to J/m3 --- !
543         DO jk = 1, nlay_i
544            WHERE( h_i_1d(1:npti)>0._wp ) e_i_1d(1:npti,jk) = e_i_1d(1:npti,jk) / (h_i_1d(1:npti) * a_i_1d(1:npti)) * nlay_i
545         END DO
546         DO jk = 1, nlay_s
547            WHERE( h_s_1d(1:npti)>0._wp ) e_s_1d(1:npti,jk) = e_s_1d(1:npti,jk) / (h_s_1d(1:npti) * a_i_1d(1:npti)) * nlay_s
548         END DO
549         !
550         !                 !---------------------!
551      CASE( 2 )            !==  from 1D to 2D  ==!
552         !                 !---------------------!
553         ! --- Change units of e_i, e_s from J/m3 to J/m2 --- !
554         DO jk = 1, nlay_i
555            e_i_1d(1:npti,jk) = e_i_1d(1:npti,jk) * h_i_1d(1:npti) * a_i_1d(1:npti) * r1_nlay_i
556         END DO
557         DO jk = 1, nlay_s
558            e_s_1d(1:npti,jk) = e_s_1d(1:npti,jk) * h_s_1d(1:npti) * a_i_1d(1:npti) * r1_nlay_s
559         END DO
560         !
561         ! Change thickness to volume (replaces routine ice_var_eqv2glo)
562         v_i_1d (1:npti) = h_i_1d (1:npti) * a_i_1d (1:npti)
563         v_s_1d (1:npti) = h_s_1d (1:npti) * a_i_1d (1:npti)
564         sv_i_1d(1:npti) = s_i_1d (1:npti) * v_i_1d (1:npti)
565         oa_i_1d(1:npti) = o_i_1d (1:npti) * a_i_1d (1:npti)
566
567         CALL tab_1d_2d( npti, nptidx(1:npti), at_i_1d(1:npti), at_i             )
568         CALL tab_1d_2d( npti, nptidx(1:npti), a_i_1d (1:npti), a_i (:,:,kl)     )
569         CALL tab_1d_2d( npti, nptidx(1:npti), h_i_1d (1:npti), h_i (:,:,kl)     )
570         CALL tab_1d_2d( npti, nptidx(1:npti), h_s_1d (1:npti), h_s (:,:,kl)     )
571         CALL tab_1d_2d( npti, nptidx(1:npti), t_su_1d(1:npti), t_su(:,:,kl)     )
572         CALL tab_1d_2d( npti, nptidx(1:npti), s_i_1d (1:npti), s_i (:,:,kl)     )
573         DO jk = 1, nlay_s
574            CALL tab_1d_2d( npti, nptidx(1:npti), t_s_1d(1:npti,jk), t_s(:,:,jk,kl)    )
575            CALL tab_1d_2d( npti, nptidx(1:npti), e_s_1d(1:npti,jk), e_s(:,:,jk,kl)    )
576         END DO
577         DO jk = 1, nlay_i
578            CALL tab_1d_2d( npti, nptidx(1:npti), t_i_1d (1:npti,jk), t_i (:,:,jk,kl)  )
579            CALL tab_1d_2d( npti, nptidx(1:npti), e_i_1d (1:npti,jk), e_i (:,:,jk,kl)  )
580            CALL tab_1d_2d( npti, nptidx(1:npti), sz_i_1d(1:npti,jk), sz_i(:,:,jk,kl)  )
581         END DO
582         !
583         CALL tab_1d_2d( npti, nptidx(1:npti), wfx_snw_sni_1d(1:npti), wfx_snw_sni )
584         CALL tab_1d_2d( npti, nptidx(1:npti), wfx_snw_sum_1d(1:npti), wfx_snw_sum )
585         CALL tab_1d_2d( npti, nptidx(1:npti), wfx_sub_1d    (1:npti), wfx_sub     )
586         CALL tab_1d_2d( npti, nptidx(1:npti), wfx_snw_sub_1d(1:npti), wfx_snw_sub )
587         CALL tab_1d_2d( npti, nptidx(1:npti), wfx_ice_sub_1d(1:npti), wfx_ice_sub )
588         CALL tab_1d_2d( npti, nptidx(1:npti), wfx_err_sub_1d(1:npti), wfx_err_sub )
589         !
590         CALL tab_1d_2d( npti, nptidx(1:npti), wfx_bog_1d (1:npti), wfx_bog        )
591         CALL tab_1d_2d( npti, nptidx(1:npti), wfx_bom_1d (1:npti), wfx_bom        )
592         CALL tab_1d_2d( npti, nptidx(1:npti), wfx_sum_1d (1:npti), wfx_sum        )
593         CALL tab_1d_2d( npti, nptidx(1:npti), wfx_sni_1d (1:npti), wfx_sni        )
594         CALL tab_1d_2d( npti, nptidx(1:npti), wfx_res_1d (1:npti), wfx_res        )
595         CALL tab_1d_2d( npti, nptidx(1:npti), wfx_spr_1d (1:npti), wfx_spr        )
596         CALL tab_1d_2d( npti, nptidx(1:npti), wfx_lam_1d (1:npti), wfx_lam        )
597         !
598         CALL tab_1d_2d( npti, nptidx(1:npti), sfx_bog_1d (1:npti), sfx_bog        )
599         CALL tab_1d_2d( npti, nptidx(1:npti), sfx_bom_1d (1:npti), sfx_bom        )
600         CALL tab_1d_2d( npti, nptidx(1:npti), sfx_sum_1d (1:npti), sfx_sum        )
601         CALL tab_1d_2d( npti, nptidx(1:npti), sfx_sni_1d (1:npti), sfx_sni        )
602         CALL tab_1d_2d( npti, nptidx(1:npti), sfx_bri_1d (1:npti), sfx_bri        )
603         CALL tab_1d_2d( npti, nptidx(1:npti), sfx_res_1d (1:npti), sfx_res        )
604         CALL tab_1d_2d( npti, nptidx(1:npti), sfx_sub_1d (1:npti), sfx_sub        )
605         CALL tab_1d_2d( npti, nptidx(1:npti), sfx_lam_1d (1:npti), sfx_lam        )
606         !
607         CALL tab_1d_2d( npti, nptidx(1:npti), hfx_thd_1d    (1:npti), hfx_thd     )
608         CALL tab_1d_2d( npti, nptidx(1:npti), hfx_spr_1d    (1:npti), hfx_spr     )
609         CALL tab_1d_2d( npti, nptidx(1:npti), hfx_sum_1d    (1:npti), hfx_sum     )
610         CALL tab_1d_2d( npti, nptidx(1:npti), hfx_bom_1d    (1:npti), hfx_bom     )
611         CALL tab_1d_2d( npti, nptidx(1:npti), hfx_bog_1d    (1:npti), hfx_bog     )
612         CALL tab_1d_2d( npti, nptidx(1:npti), hfx_dif_1d    (1:npti), hfx_dif     )
613         CALL tab_1d_2d( npti, nptidx(1:npti), hfx_opw_1d    (1:npti), hfx_opw     )
614         CALL tab_1d_2d( npti, nptidx(1:npti), hfx_snw_1d    (1:npti), hfx_snw     )
615         CALL tab_1d_2d( npti, nptidx(1:npti), hfx_sub_1d    (1:npti), hfx_sub     )
616         CALL tab_1d_2d( npti, nptidx(1:npti), hfx_res_1d    (1:npti), hfx_res     )
617         CALL tab_1d_2d( npti, nptidx(1:npti), hfx_err_dif_1d(1:npti), hfx_err_dif )
618         !
619         CALL tab_1d_2d( npti, nptidx(1:npti), qns_ice_1d    (1:npti), qns_ice    (:,:,kl) )
620         CALL tab_1d_2d( npti, nptidx(1:npti), qtr_ice_bot_1d(1:npti), qtr_ice_bot(:,:,kl) )
621         ! effective conductivity and 1st layer temperature (ln_cndflx=T)
622         CALL tab_1d_2d( npti, nptidx(1:npti), cnd_ice_1d(1:npti), cnd_ice(:,:,kl) )
623         CALL tab_1d_2d( npti, nptidx(1:npti), t1_ice_1d (1:npti), t1_ice (:,:,kl) )
624         ! Melt ponds
625         CALL tab_1d_2d( npti, nptidx(1:npti), dh_i_sum  (1:npti) , dh_i_sum_2d(:,:,kl) )
626         CALL tab_1d_2d( npti, nptidx(1:npti), dh_s_mlt  (1:npti) , dh_s_mlt_2d(:,:,kl) )
627         ! SIMIP diagnostics
628         CALL tab_1d_2d( npti, nptidx(1:npti), t_si_1d       (1:npti), t_si       (:,:,kl) )
629         CALL tab_1d_2d( npti, nptidx(1:npti), qcn_ice_bot_1d(1:npti), qcn_ice_bot(:,:,kl) )
630         CALL tab_1d_2d( npti, nptidx(1:npti), qcn_ice_top_1d(1:npti), qcn_ice_top(:,:,kl) )
631         CALL tab_1d_2d( npti, nptidx(1:npti), qml_ice_1d    (1:npti), qml_ice    (:,:,kl) )
632         ! extensive variables
633         CALL tab_1d_2d( npti, nptidx(1:npti), v_i_1d (1:npti), v_i (:,:,kl) )
634         CALL tab_1d_2d( npti, nptidx(1:npti), v_s_1d (1:npti), v_s (:,:,kl) )
635         CALL tab_1d_2d( npti, nptidx(1:npti), sv_i_1d(1:npti), sv_i(:,:,kl) )
636         CALL tab_1d_2d( npti, nptidx(1:npti), oa_i_1d(1:npti), oa_i(:,:,kl) )
637         ! check convergence of heat diffusion scheme
638         IF( ln_zdf_chkcvg ) THEN
639            CALL tab_1d_2d( npti, nptidx(1:npti), tice_cvgerr_1d(1:npti), ztice_cvgerr(:,:,kl) )
640            CALL tab_1d_2d( npti, nptidx(1:npti), tice_cvgstp_1d(1:npti), ztice_cvgstp(:,:,kl) )
641         ENDIF
642         !
643      END SELECT
644      !
645   END SUBROUTINE ice_thd_1d2d
646
647
648   SUBROUTINE ice_thd_init
649      !!-------------------------------------------------------------------
650      !!                   ***  ROUTINE ice_thd_init ***
651      !!
652      !! ** Purpose :   Physical constants and parameters associated with
653      !!                ice thermodynamics
654      !!
655      !! ** Method  :   Read the namthd namelist and check the parameters
656      !!                called at the first timestep (nit000)
657      !!
658      !! ** input   :   Namelist namthd
659      !!-------------------------------------------------------------------
660      INTEGER  ::   ios   ! Local integer output status for namelist read
661      !!
662      NAMELIST/namthd/ ln_icedH, ln_icedA, ln_icedO, ln_icedS, ln_leadhfx
663      !!-------------------------------------------------------------------
664      !
665      READ  ( numnam_ice_ref, namthd, IOSTAT = ios, ERR = 901)
666901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namthd in reference namelist' )
667      READ  ( numnam_ice_cfg, namthd, IOSTAT = ios, ERR = 902 )
668902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namthd in configuration namelist' )
669      IF(lwm) WRITE( numoni, namthd )
670      !
671      IF(lwp) THEN                          ! control print
672         WRITE(numout,*)
673         WRITE(numout,*) 'ice_thd_init: Ice Thermodynamics'
674         WRITE(numout,*) '~~~~~~~~~~~~'
675         WRITE(numout,*) '   Namelist namthd:'
676         WRITE(numout,*) '      activate ice thick change from top/bot (T) or not (F)                ln_icedH   = ', ln_icedH
677         WRITE(numout,*) '      activate lateral melting (T) or not (F)                              ln_icedA   = ', ln_icedA
678         WRITE(numout,*) '      activate ice growth in open-water (T) or not (F)                     ln_icedO   = ', ln_icedO
679         WRITE(numout,*) '      activate gravity drainage and flushing (T) or not (F)                ln_icedS   = ', ln_icedS
680         WRITE(numout,*) '      heat in the leads is used to melt sea-ice before warming the ocean   ln_leadhfx = ', ln_leadhfx
681     ENDIF
682      !
683                       CALL ice_thd_zdf_init   ! set ice heat diffusion parameters
684      IF( ln_icedA )   CALL ice_thd_da_init    ! set ice lateral melting parameters
685      IF( ln_icedO )   CALL ice_thd_do_init    ! set ice growth in open water parameters
686                       CALL ice_thd_sal_init   ! set ice salinity parameters
687                       CALL ice_thd_pnd_init   ! set melt ponds parameters
688      !
689   END SUBROUTINE ice_thd_init
690
691#else
692   !!----------------------------------------------------------------------
693   !!   Default option         Dummy module          NO  SI3 sea-ice model
694   !!----------------------------------------------------------------------
695#endif
696
697   !!======================================================================
698END MODULE icethd
Note: See TracBrowser for help on using the repository browser.