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.
icesbc.F90 in NEMO/trunk/src/ICE – NEMO

source: NEMO/trunk/src/ICE/icesbc.F90

Last change on this file was 15388, checked in by clem, 2 years ago

slightly rearrange ice thermo. No change in sette

  • Property svn:keywords set to Id
File size: 24.9 KB
Line 
1MODULE icesbc
2   !!======================================================================
3   !!                       ***  MODULE  icesbc  ***
4   !! Sea-Ice :   air-ice sbc fields
5   !!=====================================================================
6   !! History :  4.0  !  2017-08  (C. Rousset)       Original code
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   USE oce            ! ocean dynamics and tracers
14   USE dom_oce        ! ocean space and time domain
15   USE ice            ! sea-ice: variables
16   USE sbc_oce        ! Surface boundary condition: ocean fields
17   USE sbc_ice        ! Surface boundary condition: ice   fields
18   USE usrdef_sbc     ! Surface boundary condition: user defined
19   USE sbcblk         ! Surface boundary condition: bulk
20   USE sbccpl         ! Surface boundary condition: coupled interface
21   USE icealb         ! sea-ice: albedo
22   !
23   USE in_out_manager ! I/O manager
24   USE iom            ! I/O manager library
25   USE lib_mpp        ! MPP library
26   USE lib_fortran    ! fortran utilities (glob_sum + no signed zero)
27   USE lbclnk         ! lateral boundary conditions (or mpp links)
28   USE timing         ! Timing
29   USE fldread        !!GS: needed by agrif
30
31   IMPLICIT NONE
32   PRIVATE
33
34   PUBLIC ice_sbc_tau   ! called by icestp.F90
35   PUBLIC ice_sbc_flx   ! called by icestp.F90
36   PUBLIC ice_sbc_init  ! called by icestp.F90
37
38   !! * Substitutions
39#  include "do_loop_substitute.h90"
40   !!----------------------------------------------------------------------
41   !! NEMO/ICE 4.0 , NEMO Consortium (2018)
42   !! $Id$
43   !! Software governed by the CeCILL license (see ./LICENSE)
44   !!----------------------------------------------------------------------
45CONTAINS
46
47   SUBROUTINE ice_sbc_tau( kt, ksbc, utau_ice, vtau_ice )
48      !!-------------------------------------------------------------------
49      !!                  ***  ROUTINE ice_sbc_tau  ***
50      !!
51      !! ** Purpose : provide surface boundary condition for sea ice (momentum)
52      !!
53      !! ** Action  : It provides the following fields:
54      !!              utau_ice, vtau_ice : surface ice stress (U- & V-points) [N/m2]
55      !!-------------------------------------------------------------------
56      INTEGER                     , INTENT(in   ) ::   kt                   ! ocean time step
57      INTEGER                     , INTENT(in   ) ::   ksbc                 ! type of sbc flux
58      REAL(wp), DIMENSION(jpi,jpj), INTENT(  out) ::   utau_ice, vtau_ice   ! air-ice stress   [N/m2]
59      !!
60      INTEGER  ::   ji, jj                 ! dummy loop index
61      REAL(wp), DIMENSION(jpi,jpj) ::   zutau_ice, zvtau_ice
62      !!-------------------------------------------------------------------
63      !
64      IF( ln_timing )   CALL timing_start('icesbc')
65      !
66      IF( kt == nit000 .AND. lwp ) THEN
67         WRITE(numout,*)
68         WRITE(numout,*)'ice_sbc_tau: Surface boundary condition for sea ice (momentum)'
69         WRITE(numout,*)'~~~~~~~~~~~~~~~'
70      ENDIF
71      !
72      SELECT CASE( ksbc )
73         CASE( jp_usr     )   ;    CALL usrdef_sbc_ice_tau( kt )                 ! user defined formulation
74      CASE( jp_blk     )
75         CALL blk_ice_1( sf(jp_wndi)%fnow(:,:,1), sf(jp_wndj)%fnow(:,:,1),   &
76            &                                      theta_air_zt(:,:), q_air_zt(:,:),   &   ! #LB: known from "sbc_oce" module...
77            &                                      sf(jp_slp )%fnow(:,:,1), u_ice, v_ice, tm_su    ,   &   ! inputs
78            &                                      putaui = utau_ice, pvtaui = vtau_ice            )       ! outputs
79 !        CASE( jp_abl     )    utau_ice & vtau_ice are computed in ablmod
80         CASE( jp_purecpl )   ;    CALL sbc_cpl_ice_tau( utau_ice , vtau_ice )   ! Coupled      formulation
81      END SELECT
82      !
83      IF( ln_mixcpl) THEN                                                        ! Case of a mixed Bulk/Coupled formulation
84                                   CALL sbc_cpl_ice_tau( zutau_ice , zvtau_ice )
85         DO_2D( 0, 0, 0, 0 )
86            utau_ice(ji,jj) = utau_ice(ji,jj) * xcplmask(ji,jj,0) + zutau_ice(ji,jj) * ( 1. - xcplmask(ji,jj,0) )
87            vtau_ice(ji,jj) = vtau_ice(ji,jj) * xcplmask(ji,jj,0) + zvtau_ice(ji,jj) * ( 1. - xcplmask(ji,jj,0) )
88         END_2D
89         CALL lbc_lnk( 'icesbc', utau_ice, 'U', -1.0_wp, vtau_ice, 'V', -1.0_wp )
90      ENDIF
91      !
92      IF( ln_timing )   CALL timing_stop('icesbc')
93      !
94   END SUBROUTINE ice_sbc_tau
95
96
97   SUBROUTINE ice_sbc_flx( kt, ksbc )
98      !!-------------------------------------------------------------------
99      !!                  ***  ROUTINE ice_sbc_flx  ***
100      !!
101      !! ** Purpose : provide surface boundary condition for sea ice (flux)
102      !!
103      !! ** Action  : It provides the following fields used in sea ice model:
104      !!                emp_oce , emp_ice                        = E-P over ocean and sea ice                    [Kg/m2/s]
105      !!                sprecip                                  = solid precipitation                           [Kg/m2/s]
106      !!                evap_ice                                 = sublimation                                   [Kg/m2/s]
107      !!                qsr_tot , qns_tot                        = solar & non solar heat flux (total)           [W/m2]
108      !!                qsr_ice , qns_ice                        = solar & non solar heat flux over ice          [W/m2]
109      !!                dqns_ice                                 = non solar  heat sensistivity                  [W/m2]
110      !!                qemp_oce, qemp_ice, qprec_ice, qevap_ice = sensible heat (associated with evap & precip) [W/m2]
111      !!            + these fields
112      !!                qsb_ice_bot                              = sensible heat at the ice bottom               [W/m2]
113      !!                fhld, qlead                              = heat budget in the leads                      [W/m2]
114      !!            + some fields that are not used outside this module:
115      !!                qla_ice                                  = latent heat flux over ice                     [W/m2]
116      !!                dqla_ice                                 = latent heat sensistivity                      [W/m2]
117      !!                tprecip                                  = total  precipitation                          [Kg/m2/s]
118      !!                alb_ice                                  = albedo above sea ice
119      !!-------------------------------------------------------------------
120      INTEGER, INTENT(in) ::   kt     ! ocean time step
121      INTEGER, INTENT(in) ::   ksbc   ! flux formulation (user defined, bulk or Pure Coupled)
122      !!--------------------------------------------------------------------
123      !
124      IF( ln_timing )   CALL timing_start('icesbc')
125
126      IF( kt == nit000 .AND. lwp ) THEN
127         WRITE(numout,*)
128         WRITE(numout,*)'ice_sbc_flx: Surface boundary condition for sea ice (flux)'
129         WRITE(numout,*)'~~~~~~~~~~~~~~~'
130      ENDIF
131      !                     !== ice albedo ==!
132      CALL ice_alb( t_su, h_i, h_s, ln_pnd_alb, a_ip_eff, h_ip, cloud_fra, alb_ice )
133      !
134      SELECT CASE( ksbc )   !== fluxes over sea ice ==!
135      !
136      CASE( jp_usr )              !--- user defined formulation
137                                  CALL usrdef_sbc_ice_flx( kt, h_s, h_i )
138      CASE( jp_blk, jp_abl )      !--- bulk formulation & ABL formulation
139                                  CALL blk_ice_2    ( t_su, h_s, h_i, alb_ice, &
140            &                                         theta_air_zt(:,:), q_air_zt(:,:),    &   ! #LB: known from "sbc_oce" module...
141            &                                         sf(jp_slp)%fnow(:,:,1), sf(jp_qlw)%fnow(:,:,1), &
142            &                                         sf(jp_prec)%fnow(:,:,1), sf(jp_snow)%fnow(:,:,1) )
143         IF( ln_mixcpl        )   CALL sbc_cpl_ice_flx( kt, picefr=at_i_b, palbi=alb_ice, psst=sst_m, pist=t_su, phs=h_s, phi=h_i )
144         IF( nn_flxdist /= -1 )   CALL ice_flx_dist   ( t_su, alb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_flxdist )
145         !                        !    compute conduction flux and surface temperature (as in Jules surface module)
146         IF( ln_cndflx .AND. .NOT.ln_cndemulate ) &
147            &                     CALL blk_ice_qcn    ( ln_virtual_itd, t_su, t_bo, h_s, h_i )
148      CASE ( jp_purecpl )         !--- coupled formulation
149                                  CALL sbc_cpl_ice_flx( kt, picefr=at_i_b, palbi=alb_ice, psst=sst_m, pist=t_su, phs=h_s, phi=h_i )
150         IF( nn_flxdist /= -1 )   CALL ice_flx_dist   ( t_su, alb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_flxdist )
151      END SELECT
152      !                     !== some fluxes at the ice-ocean interface and in the leads
153      CALL ice_flx_other
154      !
155      IF( ln_timing )   CALL timing_stop('icesbc')
156      !
157   END SUBROUTINE ice_sbc_flx
158
159
160   SUBROUTINE ice_flx_dist( ptn_ice, palb_ice, pqns_ice, pqsr_ice, pdqn_ice, pevap_ice, pdevap_ice, k_flxdist )
161      !!-------------------------------------------------------------------
162      !!                  ***  ROUTINE ice_flx_dist  ***
163      !!
164      !! ** Purpose :   update the ice surface boundary condition by averaging
165      !!              and/or redistributing fluxes on ice categories
166      !!
167      !! ** Method  :   average then redistribute
168      !!
169      !! ** Action  :   depends on k_flxdist
170      !!                = -1  Do nothing (needs N(cat) fluxes)
171      !!                =  0  Average N(cat) fluxes then apply the average over the N(cat) ice
172      !!                =  1  Average N(cat) fluxes then redistribute over the N(cat) ice
173      !!                                                 using T-ice and albedo sensitivity
174      !!                =  2  Redistribute a single flux over categories
175      !!-------------------------------------------------------------------
176      INTEGER                   , INTENT(in   ) ::   k_flxdist  ! redistributor
177      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   ptn_ice    ! ice surface temperature
178      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   palb_ice   ! ice albedo
179      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pqns_ice   ! non solar flux
180      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pqsr_ice   ! net solar flux
181      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pdqn_ice   ! non solar flux sensitivity
182      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pevap_ice  ! sublimation
183      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pdevap_ice ! sublimation sensitivity
184      !
185      INTEGER  ::   jl      ! dummy loop index
186      !
187      REAL(wp), DIMENSION(jpi,jpj) ::   z1_at_i   ! inverse of concentration
188      !
189      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   z_qsr_m   ! Mean solar heat flux over all categories
190      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   z_qns_m   ! Mean non solar heat flux over all categories
191      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   z_evap_m  ! Mean sublimation over all categories
192      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   z_dqn_m   ! Mean d(qns)/dT over all categories
193      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   z_devap_m ! Mean d(evap)/dT over all categories
194      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   zalb_m    ! Mean albedo over all categories
195      REAL(wp), ALLOCATABLE, DIMENSION(:,:) ::   ztem_m    ! Mean temperature over all categories
196      !!----------------------------------------------------------------------
197      !
198      WHERE ( at_i (:,:) > 0._wp )   ; z1_at_i(:,:) = 1._wp / at_i (:,:)
199      ELSEWHERE                      ; z1_at_i(:,:) = 0._wp
200      END WHERE
201
202      SELECT CASE( k_flxdist )       !==  averaged on all ice categories  ==!
203      !
204      CASE( 0 , 1 )
205         !
206         ALLOCATE( z_qns_m(jpi,jpj), z_qsr_m(jpi,jpj), z_dqn_m(jpi,jpj), z_evap_m(jpi,jpj), z_devap_m(jpi,jpj) )
207         !
208         z_qns_m  (:,:) = SUM( a_i(:,:,:) * pqns_ice  (:,:,:) , dim=3 ) * z1_at_i(:,:)
209         z_qsr_m  (:,:) = SUM( a_i(:,:,:) * pqsr_ice  (:,:,:) , dim=3 ) * z1_at_i(:,:)
210         z_dqn_m  (:,:) = SUM( a_i(:,:,:) * pdqn_ice  (:,:,:) , dim=3 ) * z1_at_i(:,:)
211         z_evap_m (:,:) = SUM( a_i(:,:,:) * pevap_ice (:,:,:) , dim=3 ) * z1_at_i(:,:)
212         z_devap_m(:,:) = SUM( a_i(:,:,:) * pdevap_ice(:,:,:) , dim=3 ) * z1_at_i(:,:)
213         DO jl = 1, jpl
214            pqns_ice  (:,:,jl) = z_qns_m (:,:)
215            pqsr_ice  (:,:,jl) = z_qsr_m (:,:)
216            pdqn_ice  (:,:,jl) = z_dqn_m  (:,:)
217            pevap_ice (:,:,jl) = z_evap_m(:,:)
218            pdevap_ice(:,:,jl) = z_devap_m(:,:)
219         END DO
220         !
221         DEALLOCATE( z_qns_m, z_qsr_m, z_dqn_m, z_evap_m, z_devap_m )
222         !
223      END SELECT
224      !
225      SELECT CASE( k_flxdist )       !==  redistribution on all ice categories  ==!
226      !
227      CASE( 1 , 2 )
228         !
229         ALLOCATE( zalb_m(jpi,jpj), ztem_m(jpi,jpj) )
230         !
231         zalb_m(:,:) = SUM( a_i(:,:,:) * palb_ice(:,:,:) , dim=3 ) * z1_at_i(:,:)
232         ztem_m(:,:) = SUM( a_i(:,:,:) * ptn_ice (:,:,:) , dim=3 ) * z1_at_i(:,:)
233         DO jl = 1, jpl
234            pqns_ice (:,:,jl) = pqns_ice (:,:,jl) + pdqn_ice  (:,:,jl) * ( ptn_ice(:,:,jl) - ztem_m(:,:) )
235            pevap_ice(:,:,jl) = pevap_ice(:,:,jl) + pdevap_ice(:,:,jl) * ( ptn_ice(:,:,jl) - ztem_m(:,:) )
236            pqsr_ice (:,:,jl) = pqsr_ice (:,:,jl) * ( 1._wp - palb_ice(:,:,jl) ) / ( 1._wp - zalb_m(:,:) )
237         END DO
238         !
239         DEALLOCATE( zalb_m, ztem_m )
240         !
241      END SELECT
242      !
243   END SUBROUTINE ice_flx_dist
244
245
246   SUBROUTINE ice_flx_other
247      !!-----------------------------------------------------------------------
248      !!                   ***  ROUTINE ice_flx_other ***
249      !!
250      !! ** Purpose :   prepare necessary fields for thermo calculations
251      !!
252      !! ** Inputs  :   u_ice, v_ice, ssu_m, ssv_m, utau, vtau
253      !!                frq_m, qsr_oce, qns_oce, qemp_oce, e3t_m, sst_m
254      !! ** Outputs :   qsb_ice_bot, fhld, qlead
255      !!-----------------------------------------------------------------------
256      INTEGER  ::   ji, jj             ! dummy loop indices
257      REAL(wp) ::   zfric_u, zqld, zqfr, zqfr_neg, zqfr_pos, zu_io, zv_io, zu_iom1, zv_iom1
258      REAL(wp), PARAMETER ::   zfric_umin = 0._wp       ! lower bound for the friction velocity (cice value=5.e-04)
259      REAL(wp), PARAMETER ::   zch        = 0.0057_wp   ! heat transfer coefficient
260      REAL(wp), DIMENSION(jpi,jpj) ::  zfric, zvel      ! ice-ocean velocity (m/s) and frictional velocity (m2/s2)
261      !!-----------------------------------------------------------------------
262      !
263      ! computation of friction velocity at T points
264      IF( ln_icedyn ) THEN
265         DO_2D( 0, 0, 0, 0 )
266            zu_io   = u_ice(ji  ,jj  ) - ssu_m(ji  ,jj  )
267            zu_iom1 = u_ice(ji-1,jj  ) - ssu_m(ji-1,jj  )
268            zv_io   = v_ice(ji  ,jj  ) - ssv_m(ji  ,jj  )
269            zv_iom1 = v_ice(ji  ,jj-1) - ssv_m(ji  ,jj-1)
270            !
271            zfric(ji,jj) = rn_cio * ( 0.5_wp * ( zu_io*zu_io + zu_iom1*zu_iom1 + zv_io*zv_io + zv_iom1*zv_iom1 ) ) * tmask(ji,jj,1)
272            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) ) + &
273               &                          ( v_ice(ji  ,jj-1) + v_ice(ji,jj) ) * ( v_ice(ji  ,jj-1) + v_ice(ji,jj) ) )
274         END_2D
275      ELSE      !  if no ice dynamics => transfer directly the atmospheric stress to the ocean
276         DO_2D( 0, 0, 0, 0 )
277            zfric(ji,jj) = r1_rho0 * SQRT( 0.5_wp *  &
278               &                         (  utau(ji,jj) * utau(ji,jj) + utau(ji-1,jj) * utau(ji-1,jj)   &
279               &                          + vtau(ji,jj) * vtau(ji,jj) + vtau(ji,jj-1) * vtau(ji,jj-1) ) ) * tmask(ji,jj,1)
280            zvel(ji,jj) = 0._wp
281         END_2D
282      ENDIF
283      CALL lbc_lnk( 'icesbc', zfric, 'T',  1.0_wp, zvel, 'T', 1.0_wp )
284      !
285      !--------------------------------------------------------------------!
286      ! Partial computation of forcing for the thermodynamic sea ice model
287      !--------------------------------------------------------------------!
288      DO_2D( nn_hls, nn_hls, nn_hls, nn_hls )   ! needed for qlead
289         rswitch  = tmask(ji,jj,1) * MAX( 0._wp , SIGN( 1._wp , at_i(ji,jj) - epsi10 ) ) ! 0 if no ice
290         !
291         ! --- Energy received in the lead from atm-oce exchanges, zqld is defined everywhere (J.m-2) --- !
292         zqld =  tmask(ji,jj,1) * rDt_ice *  &
293            &    ( ( 1._wp - at_i_b(ji,jj) ) * qsr_oce(ji,jj) * frq_m(ji,jj) +  &
294            &      ( 1._wp - at_i_b(ji,jj) ) * qns_oce(ji,jj) + qemp_oce(ji,jj) )
295
296         ! --- Energy needed to bring ocean surface layer until its freezing, zqfr is defined everywhere (J.m-2) --- !
297         !     (mostly<0 but >0 if supercooling)
298         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)
299         zqfr_neg = MIN( zqfr , 0._wp )                                                                    ! only < 0
300         zqfr_pos = MAX( zqfr , 0._wp )                                                                    ! only > 0
301
302         ! --- Sensible ocean-to-ice heat flux (W/m2) --- !
303         !     (mostly>0 but <0 if supercooling)
304         zfric_u            = MAX( SQRT( zfric(ji,jj) ), zfric_umin )
305         qsb_ice_bot(ji,jj) = rswitch * rho0 * rcp * zch * zfric_u * ( ( sst_m(ji,jj) + rt0 ) - t_bo(ji,jj) )
306
307         ! upper bound for qsb_ice_bot: the heat retrieved from the ocean must be smaller than the heat necessary to reach
308         !                              the freezing point, so that we do not have SST < T_freeze
309         !                              This implies: qsb_ice_bot(ji,jj) * at_i(ji,jj) * rtdice <= - zqfr_neg
310         !                              The following formulation is ok for both normal conditions and supercooling
311         qsb_ice_bot(ji,jj) = rswitch * MIN( qsb_ice_bot(ji,jj), - zqfr_neg * r1_Dt_ice / MAX( at_i(ji,jj), epsi10 ) )
312
313         ! If conditions are always supercooled (such as at the mouth of ice-shelves), then ice grows continuously
314         ! ==> stop ice formation by artificially setting up the turbulent fluxes to 0 when volume > 20m (arbitrary)
315         IF( ( t_bo(ji,jj) - ( sst_m(ji,jj) + rt0 ) ) > 0._wp .AND. vt_i(ji,jj) >= 20._wp ) THEN
316            zqfr               = 0._wp
317            zqfr_pos           = 0._wp
318            qsb_ice_bot(ji,jj) = 0._wp
319         ENDIF
320         !
321         ! --- Energy Budget of the leads (qlead, J.m-2) --- !
322         !     qlead is the energy received from the atm. in the leads.
323         !     If warming (zqld >= 0), then the energy in the leads is used to melt ice (bottom melting) => fhld  (W/m2)
324         !     If cooling (zqld <  0), then the energy in the leads is used to grow ice in open water    => qlead (J.m-2)
325         IF( zqld >= 0._wp .AND. at_i(ji,jj) > 0._wp ) THEN
326            ! upper bound for fhld: fhld should be equal to zqld
327            !                        but we have to make sure that this heat will not make the sst drop below the freezing point
328            !                        so the max heat that can be pulled out of the ocean is zqld - qsb - zqfr_pos
329            !                        The following formulation is ok for both normal conditions and supercooling
330            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
331               &                                 - qsb_ice_bot(ji,jj) )
332            qlead(ji,jj) = 0._wp
333         ELSE
334            fhld (ji,jj) = 0._wp
335            ! upper bound for qlead: qlead should be equal to zqld
336            !                        but before using this heat for ice formation, we suppose that the ocean cools down till the freezing point.
337            !                        The energy for this cooling down is zqfr. Also some heat will be removed from the ocean from turbulent fluxes (qsb)
338            !                        and freezing point is reached if zqfr = zqld - qsb*a/dt
339            !                        so the max heat that can be pulled out of the ocean is zqld - qsb - zqfr
340            !                        The following formulation is ok for both normal conditions and supercooling
341            qlead(ji,jj) = MIN( 0._wp , zqld - ( qsb_ice_bot(ji,jj) * at_i(ji,jj) * rDt_ice ) - zqfr )
342         ENDIF
343         !
344         ! If ice is landfast and ice concentration reaches its max
345         ! => stop ice formation in open water
346         IF(  zvel(ji,jj) <= 5.e-04_wp .AND. at_i(ji,jj) >= rn_amax_2d(ji,jj)-epsi06 )   qlead(ji,jj) = 0._wp
347         !
348         ! If the grid cell is almost fully covered by ice (no leads)
349         ! => stop ice formation in open water
350         IF( at_i(ji,jj) >= (1._wp - epsi10) )   qlead(ji,jj) = 0._wp
351         !
352         ! If ln_leadhfx is false
353         ! => do not use energy of the leads to melt sea-ice
354         IF( .NOT.ln_leadhfx )   fhld(ji,jj) = 0._wp
355         !
356      END_2D
357
358      ! In case we bypass open-water ice formation
359      IF( .NOT. ln_icedO )  qlead(:,:) = 0._wp
360      ! In case we bypass growing/melting from top and bottom
361      IF( .NOT. ln_icedH ) THEN
362         qsb_ice_bot(:,:) = 0._wp
363         fhld       (:,:) = 0._wp
364      ENDIF
365     
366   END SUBROUTINE ice_flx_other
367   
368   
369   SUBROUTINE ice_sbc_init
370      !!-------------------------------------------------------------------
371      !!                  ***  ROUTINE ice_sbc_init  ***
372      !!
373      !! ** Purpose :   Physical constants and parameters linked to the ice dynamics
374      !!
375      !! ** Method  :   Read the namsbc namelist and check the ice-dynamic
376      !!              parameter values called at the first timestep (nit000)
377      !!
378      !! ** input   :   Namelist namsbc
379      !!-------------------------------------------------------------------
380      INTEGER ::   ios, ioptio   ! Local integer
381      !!
382      NAMELIST/namsbc/ rn_cio, nn_snwfra, rn_snwblow, nn_flxdist, ln_cndflx, ln_cndemulate, nn_qtrice
383      !!-------------------------------------------------------------------
384      !
385      READ  ( numnam_ice_ref, namsbc, IOSTAT = ios, ERR = 901)
386901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namsbc in reference namelist' )
387      READ  ( numnam_ice_cfg, namsbc, IOSTAT = ios, ERR = 902 )
388902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namsbc in configuration namelist' )
389      IF(lwm) WRITE( numoni, namsbc )
390      !
391      IF(lwp) THEN                     ! control print
392         WRITE(numout,*)
393         WRITE(numout,*) 'ice_sbc_init: ice parameters for ice dynamics '
394         WRITE(numout,*) '~~~~~~~~~~~~~~~~'
395         WRITE(numout,*) '   Namelist namsbc:'
396         WRITE(numout,*) '      drag coefficient for oceanic stress                       rn_cio        = ', rn_cio
397         WRITE(numout,*) '      fraction of ice covered by snow (options 0,1,2)           nn_snwfra     = ', nn_snwfra
398         WRITE(numout,*) '      coefficient for ice-lead partition of snowfall            rn_snwblow    = ', rn_snwblow
399         WRITE(numout,*) '      Multicategory heat flux formulation                       nn_flxdist    = ', nn_flxdist
400         WRITE(numout,*) '      Use conduction flux as surface condition                  ln_cndflx     = ', ln_cndflx
401         WRITE(numout,*) '         emulate conduction flux                                ln_cndemulate = ', ln_cndemulate
402         WRITE(numout,*) '      solar flux transmitted thru the surface scattering layer  nn_qtrice     = ', nn_qtrice
403         WRITE(numout,*) '         = 0  Grenfell and Maykut 1977'
404         WRITE(numout,*) '         = 1  Lebrun 2019'
405      ENDIF
406      !
407      IF(lwp) WRITE(numout,*)
408      SELECT CASE( nn_flxdist )         ! SI3 Multi-category heat flux formulation
409      CASE( -1  )
410         IF(lwp) WRITE(numout,*) '   SI3: use per-category fluxes (nn_flxdist = -1) '
411      CASE(  0  )
412         IF(lwp) WRITE(numout,*) '   SI3: use average per-category fluxes (nn_flxdist = 0) '
413      CASE(  1  )
414         IF(lwp) WRITE(numout,*) '   SI3: use average then redistribute per-category fluxes (nn_flxdist = 1) '
415         IF( ln_cpl )         CALL ctl_stop( 'ice_thd_init: the chosen nn_flxdist for SI3 in coupled mode must be /=1' )
416      CASE(  2  )
417         IF(lwp) WRITE(numout,*) '   SI3: Redistribute a single flux over categories (nn_flxdist = 2) '
418         IF( .NOT. ln_cpl )   CALL ctl_stop( 'ice_thd_init: the chosen nn_flxdist for SI3 in forced mode must be /=2' )
419      CASE DEFAULT
420         CALL ctl_stop( 'ice_thd_init: SI3 option, nn_flxdist, should be between -1 and 2' )
421      END SELECT
422      !
423   END SUBROUTINE ice_sbc_init
424
425#else
426   !!----------------------------------------------------------------------
427   !!   Default option :         Empty module         NO SI3 sea-ice model
428   !!----------------------------------------------------------------------
429#endif
430
431   !!======================================================================
432END MODULE icesbc
Note: See TracBrowser for help on using the repository browser.