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.
sbcice_lim.F90 in branches/2014/dev_r4621_NOC4_BDY_VERT_INTERP/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2014/dev_r4621_NOC4_BDY_VERT_INTERP/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90 @ 5038

Last change on this file since 5038 was 5038, checked in by jamesharle, 9 years ago

Merging branch with HEAD of the trunk

  • Property svn:keywords set to Id
File size: 45.5 KB
Line 
1MODULE sbcice_lim
2   !!======================================================================
3   !!                       ***  MODULE  sbcice_lim  ***
4   !! Surface module :  update the ocean surface boundary condition over ice
5   !!       &           covered area using LIM sea-ice model
6   !! Sea-Ice model  :  LIM-3 Sea ice model time-stepping
7   !!=====================================================================
8   !! History :  2.0  ! 2006-12  (M. Vancoppenolle) Original code
9   !!            3.0  ! 2008-02  (C. Talandier)  Surface module from icestp.F90
10   !!             -   ! 2008-04  (G. Madec)  sltyle and lim_ctl routine
11   !!            3.3  ! 2010-11  (G. Madec) ice-ocean stress always computed at each ocean time-step
12   !!            3.4  ! 2011-01  (A Porter)  dynamical allocation
13   !!             -   ! 2012-10  (C. Rousset)  add lim_diahsb
14   !!            3.6  ! 2014-07  (M. Vancoppenolle, G. Madec, O. Marti) revise coupled interface
15   !!----------------------------------------------------------------------
16#if defined key_lim3
17   !!----------------------------------------------------------------------
18   !!   'key_lim3' :                                  LIM 3.0 sea-ice model
19   !!----------------------------------------------------------------------
20   !!   sbc_ice_lim  : sea-ice model time-stepping and update ocean sbc over ice-covered area
21   !!   lim_ctl       : alerts in case of ice model crash
22   !!   lim_prt_state : ice control print at a given grid point
23   !!----------------------------------------------------------------------
24   USE oce             ! ocean dynamics and tracers
25   USE dom_oce         ! ocean space and time domain
26   USE par_ice         ! sea-ice parameters
27   USE ice             ! LIM-3: ice variables
28   USE iceini          ! LIM-3: ice initialisation
29   USE dom_ice         ! LIM-3: ice domain
30
31   USE sbc_oce         ! Surface boundary condition: ocean fields
32   USE sbc_ice         ! Surface boundary condition: ice   fields
33   USE sbcblk_core     ! Surface boundary condition: CORE bulk
34   USE sbcblk_clio     ! Surface boundary condition: CLIO bulk
35   USE sbccpl          ! Surface boundary condition: coupled interface
36   USE albedo          ! ocean & ice albedo
37
38   USE phycst          ! Define parameters for the routines
39   USE eosbn2          ! equation of state
40   USE limdyn          ! Ice dynamics
41   USE limtrp          ! Ice transport
42   USE limthd          ! Ice thermodynamics
43   USE limitd_th       ! Thermodynamics on ice thickness distribution
44   USE limitd_me       ! Mechanics on ice thickness distribution
45   USE limsbc          ! sea surface boundary condition
46   USE limdiahsb       ! Ice budget diagnostics
47   USE limwri          ! Ice outputs
48   USE limrst          ! Ice restarts
49   USE limupdate1       ! update of global variables
50   USE limupdate2       ! update of global variables
51   USE limvar          ! Ice variables switch
52
53   USE c1d             ! 1D vertical configuration
54   USE lbclnk          ! lateral boundary condition - MPP link
55   USE lib_mpp         ! MPP library
56   USE wrk_nemo        ! work arrays
57   USE timing          ! Timing
58   USE iom             ! I/O manager library
59   USE in_out_manager  ! I/O manager
60   USE prtctl          ! Print control
61   USE lib_fortran     !
62   USE cpl_oasis3, ONLY : lk_cpl
63
64#if defined key_bdy 
65   USE bdyice_lim       ! unstructured open boundary data  (bdy_ice_lim routine)
66#endif
67
68   IMPLICIT NONE
69   PRIVATE
70
71   PUBLIC sbc_ice_lim  ! routine called by sbcmod.F90
72   PUBLIC lim_prt_state
73   
74   !! * Substitutions
75#  include "domzgr_substitute.h90"
76#  include "vectopt_loop_substitute.h90"
77   !!----------------------------------------------------------------------
78   !! NEMO/OPA 4.0 , UCL NEMO Consortium (2011)
79   !! $Id$
80   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
81   !!----------------------------------------------------------------------
82CONTAINS
83
84   !!======================================================================
85
86   SUBROUTINE sbc_ice_lim( kt, kblk )
87      !!---------------------------------------------------------------------
88      !!                  ***  ROUTINE sbc_ice_lim  ***
89      !!                   
90      !! ** Purpose :   update the ocean surface boundary condition via the
91      !!                Louvain la Neuve Sea Ice Model time stepping
92      !!
93      !! ** Method  :   ice model time stepping
94      !!              - call the ice dynamics routine
95      !!              - call the ice advection/diffusion routine
96      !!              - call the ice thermodynamics routine
97      !!              - call the routine that computes mass and
98      !!                heat fluxes at the ice/ocean interface
99      !!              - save the outputs
100      !!              - save the outputs for restart when necessary
101      !!
102      !! ** Action  : - time evolution of the LIM sea-ice model
103      !!              - update all sbc variables below sea-ice:
104      !!                utau, vtau, taum, wndm, qns , qsr, emp , sfx
105      !!---------------------------------------------------------------------
106      INTEGER, INTENT(in) ::   kt      ! ocean time step
107      INTEGER, INTENT(in) ::   kblk    ! type of bulk (=3 CLIO, =4 CORE, =5 COUPLED)
108      !!
109      INTEGER  ::   ji, jj, jl, jk      ! dummy loop index
110      REAL(wp) ::   zcoef   ! local scalar
111      REAL(wp), POINTER, DIMENSION(:,:,:)   ::   zalb_os, zalb_cs  ! ice albedo under overcast/clear sky
112      REAL(wp), POINTER, DIMENSION(:,:,:)   ::   zalb_ice          ! mean ice albedo (for coupled)
113      !!----------------------------------------------------------------------
114
115      IF( nn_timing == 1 )  CALL timing_start('sbc_ice_lim')
116
117      IF( kt == nit000 ) THEN
118         IF(lwp) WRITE(numout,*)
119         IF(lwp) WRITE(numout,*) 'sbc_ice_lim : update ocean surface boudary condition' 
120         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   via Louvain la Neuve Ice Model (LIM-3) time stepping'
121         !
122         CALL ice_init
123         !
124         IF( ln_nicep ) THEN      ! control print at a given point
125            jiindx = 15    ;   jjindx =  44
126            IF(lwp) WRITE(numout,*) ' The debugging point is : jiindx : ',jiindx, ' jjindx : ',jjindx
127         ENDIF
128      ENDIF
129
130      !                                        !----------------------!
131      IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN     !  Ice time-step only  !
132         !                                     !----------------------!
133         !                                           !  Bulk Formulae !
134         !                                           !----------------!
135         !
136         u_oce(:,:) = ssu_m(:,:) * umask(:,:,1)                     ! mean surface ocean current at ice velocity point
137         v_oce(:,:) = ssv_m(:,:) * vmask(:,:,1)                    ! (C-grid dynamics :  U- & V-points as the ocean)
138         !
139         t_bo(:,:) = ( eos_fzp( sss_m ) +  rt0 ) * tmask(:,:,1) + rt0 * ( 1. - tmask(:,:,1) )  ! masked sea surface freezing temperature [Kelvin]
140         !                                                                                  ! (set to rt0 over land)
141         !                                           ! Ice albedo
142         CALL wrk_alloc( jpi,jpj,jpl, zalb_os, zalb_cs, zalb_ice )     
143
144         CALL albedo_ice( t_su, ht_i, ht_s, zalb_cs, zalb_os )  ! cloud-sky and overcast-sky ice albedos
145
146         SELECT CASE( kblk )
147         CASE( jp_core , jp_cpl )   ! CORE and COUPLED bulk formulations
148
149            ! albedo depends on cloud fraction because of non-linear spectral effects
150            zalb_ice(:,:,:) = ( 1. - cldf_ice ) * zalb_cs(:,:,:) + cldf_ice * zalb_os(:,:,:)
151            ! In CLIO the cloud fraction is read in the climatology and the all-sky albedo
152            ! (zalb_ice) is computed within the bulk routine
153           
154         END SELECT
155         
156         !                                           ! Mask sea ice surface temperature
157         DO jl = 1, jpl
158            t_su(:,:,jl) = t_su(:,:,jl) +  rt0 * ( 1. - tmask(:,:,1) )
159         END DO
160     
161         ! Bulk formulae  - provides the following fields:
162         ! utau_ice, vtau_ice : surface ice stress                     (U- & V-points)   [N/m2]
163         ! qsr_ice , qns_ice  : solar & non solar heat flux over ice   (T-point)         [W/m2]
164         ! qla_ice            : latent heat flux over ice              (T-point)         [W/m2]
165         ! dqns_ice, dqla_ice : non solar & latent heat sensistivity   (T-point)         [W/m2]
166         ! tprecip , sprecip  : total & solid precipitation            (T-point)         [Kg/m2/s]
167         ! fr1_i0  , fr2_i0   : 1sr & 2nd fraction of qsr penetration in ice             [%]
168         !
169         SELECT CASE( kblk )
170         CASE( jp_clio )                                       ! CLIO bulk formulation
171            CALL blk_ice_clio( t_su , zalb_cs    , zalb_os    , zalb_ice  ,               &
172               &                      utau_ice   , vtau_ice   , qns_ice   , qsr_ice   ,   &
173               &                      qla_ice    , dqns_ice   , dqla_ice  ,               &
174               &                      tprecip    , sprecip    ,                           &
175               &                      fr1_i0     , fr2_i0     , cp_ice_msh, jpl  )
176            !         
177            IF( nn_limflx /= 2 )   CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice ,   &
178               &                                           dqns_ice, qla_ice, dqla_ice, nn_limflx )
179
180         CASE( jp_core )                                       ! CORE bulk formulation
181            CALL blk_ice_core( t_su , u_ice     , v_ice     , zalb_ice   ,               &
182               &                      utau_ice  , vtau_ice  , qns_ice    , qsr_ice   ,   &
183               &                      qla_ice   , dqns_ice  , dqla_ice   ,               &
184               &                      tprecip   , sprecip   ,                            &
185               &                      fr1_i0    , fr2_i0    , cp_ice_msh, jpl  )
186               !
187            IF( nn_limflx /= 2 )   CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice ,   &
188               &                                           dqns_ice, qla_ice, dqla_ice, nn_limflx )
189            !
190         CASE ( jp_cpl )
191           
192            CALL sbc_cpl_ice_tau( utau_ice , vtau_ice )
193
194            ! MV -> seb
195!           CALL sbc_cpl_ice_flx( p_frld=ato_i, palbi=zalb_ice, psst=sst_m, pist=t_su    )
196
197!           IF( nn_limflx == 2 )   CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice ,   &
198!              &                                           dqns_ice, qla_ice, dqla_ice, nn_limflx )
199!           ! Latent heat flux is forced to 0 in coupled :
200!           !  it is included in qns (non-solar heat flux)
201!           qla_ice  (:,:,:) = 0._wp
202!           dqla_ice (:,:,:) = 0._wp
203            ! END MV -> seb
204            !
205         END SELECT
206         
207         !                                           !----------------------!
208         !                                           ! LIM-3  time-stepping !
209         !                                           !----------------------!
210         !
211         numit = numit + nn_fsbc                     ! Ice model time step
212         !
213         !                                           ! Store previous ice values
214         a_i_b  (:,:,:)   = a_i  (:,:,:)     ! ice area
215         e_i_b  (:,:,:,:) = e_i  (:,:,:,:)   ! ice thermal energy
216         v_i_b  (:,:,:)   = v_i  (:,:,:)     ! ice volume
217         v_s_b  (:,:,:)   = v_s  (:,:,:)     ! snow volume
218         e_s_b  (:,:,:,:) = e_s  (:,:,:,:)   ! snow thermal energy
219         smv_i_b(:,:,:)   = smv_i(:,:,:)     ! salt content
220         oa_i_b (:,:,:)   = oa_i (:,:,:)     ! areal age content
221         u_ice_b(:,:)     = u_ice(:,:)
222         v_ice_b(:,:)     = v_ice(:,:)
223
224         ! salt, heat and mass fluxes
225         sfx    (:,:) = 0._wp   ;
226         sfx_bri(:,:) = 0._wp   ; 
227         sfx_sni(:,:) = 0._wp   ;   sfx_opw(:,:) = 0._wp
228         sfx_bog(:,:) = 0._wp   ;   sfx_dyn(:,:) = 0._wp
229         sfx_bom(:,:) = 0._wp   ;   sfx_sum(:,:) = 0._wp
230         sfx_res(:,:) = 0._wp
231
232         wfx_snw(:,:) = 0._wp   ;   wfx_ice(:,:) = 0._wp
233         wfx_sni(:,:) = 0._wp   ;   wfx_opw(:,:) = 0._wp
234         wfx_bog(:,:) = 0._wp   ;   wfx_dyn(:,:) = 0._wp
235         wfx_bom(:,:) = 0._wp   ;   wfx_sum(:,:) = 0._wp
236         wfx_res(:,:) = 0._wp   ;   wfx_sub(:,:) = 0._wp
237         wfx_spr(:,:) = 0._wp   ;   
238
239         hfx_in (:,:) = 0._wp   ;   hfx_out(:,:) = 0._wp
240         hfx_thd(:,:) = 0._wp   ;   
241         hfx_snw(:,:) = 0._wp   ;   hfx_opw(:,:) = 0._wp
242         hfx_bog(:,:) = 0._wp   ;   hfx_dyn(:,:) = 0._wp
243         hfx_bom(:,:) = 0._wp   ;   hfx_sum(:,:) = 0._wp
244         hfx_res(:,:) = 0._wp   ;   hfx_sub(:,:) = 0._wp
245         hfx_spr(:,:) = 0._wp   ;   hfx_dif(:,:) = 0._wp 
246         hfx_err(:,:) = 0._wp   ;   hfx_err_rem(:,:) = 0._wp
247
248                          CALL lim_rst_opn( kt )     ! Open Ice restart file
249         !
250         IF( ln_nicep )   CALL lim_prt_state( kt, jiindx, jjindx, 1, ' - Beginning the time step - ' )   ! control print
251         ! ----------------------------------------------
252         ! ice dynamics and transport (except in 1D case)
253         ! ----------------------------------------------
254         IF( .NOT. lk_c1d ) THEN
255                          CALL lim_dyn( kt )              ! Ice dynamics    ( rheology/dynamics )
256                          CALL lim_trp( kt )              ! Ice transport   ( Advection/diffusion )
257                          CALL lim_var_glo2eqv            ! equivalent variables, requested for rafting
258         IF( ln_nicep )   CALL lim_prt_state( kt, jiindx, jjindx,-1, ' - ice dyn & trp - ' )   ! control print
259                          CALL lim_itd_me                 ! Mechanical redistribution ! (ridging/rafting)
260                          CALL lim_var_agg( 1 ) 
261#if defined key_bdy
262                          ! bdy ice thermo
263                          CALL lim_var_glo2eqv            ! equivalent variables
264                          CALL bdy_ice_lim( kt )
265                          CALL lim_itd_me_zapsmall
266                          CALL lim_var_agg(1)
267         IF( ln_nicep )   CALL lim_prt_state( kt, jiindx, jjindx, 1, ' - ice thermo bdy - ' )   ! control print
268#endif
269                          CALL lim_update1
270         ENDIF
271!                         !- Change old values for new values
272                          u_ice_b(:,:)     = u_ice(:,:)
273                          v_ice_b(:,:)     = v_ice(:,:)
274                          a_i_b  (:,:,:)   = a_i  (:,:,:)
275                          v_s_b  (:,:,:)   = v_s  (:,:,:)
276                          v_i_b  (:,:,:)   = v_i  (:,:,:)
277                          e_s_b  (:,:,:,:) = e_s  (:,:,:,:)
278                          e_i_b  (:,:,:,:) = e_i  (:,:,:,:)
279                          oa_i_b (:,:,:)   = oa_i (:,:,:)
280                          smv_i_b(:,:,:)   = smv_i(:,:,:)
281 
282         ! ----------------------------------------------
283         ! ice thermodynamic
284         ! ----------------------------------------------
285                          CALL lim_var_glo2eqv            ! equivalent variables
286                          CALL lim_var_agg(1)             ! aggregate ice categories
287                          ! previous lead fraction and ice volume for flux calculations
288                          pfrld(:,:)   = 1._wp - at_i(:,:)
289                          phicif(:,:)  = vt_i(:,:)
290
291                          ! MV -> seb
292                          SELECT CASE( kblk )
293                             CASE ( jp_cpl )
294                             CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=zalb_ice, psst=sst_m, pist=t_su    )
295                             IF( nn_limflx == 2 )   CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice ,   &
296                          &                                           dqns_ice, qla_ice, dqla_ice, nn_limflx )
297                           ! Latent heat flux is forced to 0 in coupled :
298                           !  it is included in qns (non-solar heat flux)
299                             qla_ice  (:,:,:) = 0._wp
300                             dqla_ice (:,:,:) = 0._wp
301                          END SELECT
302                          ! END MV -> seb
303                          !
304                          CALL lim_var_bv                 ! bulk brine volume (diag)
305                          CALL lim_thd( kt )              ! Ice thermodynamics
306                          zcoef = rdt_ice /rday           !  Ice natural aging
307                          oa_i(:,:,:) = oa_i(:,:,:) + a_i(:,:,:) * zcoef
308         IF( ln_nicep )   CALL lim_prt_state( kt, jiindx, jjindx, 1, ' - ice thermodyn. - ' )   ! control print
309                          CALL lim_itd_th( kt )           !  Remap ice categories, lateral accretion  !
310                          CALL lim_var_agg( 1 )           ! requested by limupdate
311                          CALL lim_update2                ! Global variables update
312
313                          CALL lim_var_glo2eqv            ! equivalent variables (outputs)
314                          CALL lim_var_agg(2)             ! aggregate ice thickness categories
315         IF( ln_nicep )   CALL lim_prt_state( kt, jiindx, jjindx, 2, ' - Final state - ' )   ! control print
316         !
317                          CALL lim_sbc_flx( kt )     ! Update surface ocean mass, heat and salt fluxes
318         !
319         IF( ln_nicep )   CALL lim_prt_state( kt, jiindx, jjindx, 3, ' - Final state lim_sbc - ' )   ! control print
320         !
321         !                                           ! Diagnostics and outputs
322         IF (ln_limdiaout) CALL lim_diahsb
323
324                          CALL lim_wri( 1  )              ! Ice outputs
325
326         IF( kt == nit000 .AND. ln_rstart )   &
327            &             CALL iom_close( numrir )        ! clem: close input ice restart file
328         !
329         IF( lrst_ice )   CALL lim_rst_write( kt )        ! Ice restart file
330                          CALL lim_var_glo2eqv            ! ???
331         !
332         IF( ln_nicep )   CALL lim_ctl( kt )              ! alerts in case of model crash
333         !
334         CALL wrk_dealloc( jpi,jpj,jpl, zalb_os, zalb_cs, zalb_ice )
335         !
336      ENDIF                                    ! End sea-ice time step only
337
338      !                                        !--------------------------!
339      !                                        !  at all ocean time step  !
340      !                                        !--------------------------!
341      !                                               
342      !                                              ! Update surface ocean stresses (only in ice-dynamic case)
343      !                                                   ! otherwise the atm.-ocean stresses are used everywhere
344      IF( ln_limdyn )     CALL lim_sbc_tau( kt, ub(:,:,1), vb(:,:,1) )  ! using before instantaneous surf. currents
345!!gm   remark, the ocean-ice stress is not saved in ice diag call above .....  find a solution!!!
346
347      !
348      IF( nn_timing == 1 )  CALL timing_stop('sbc_ice_lim')
349      !
350   END SUBROUTINE sbc_ice_lim
351   
352   
353      SUBROUTINE ice_lim_flx( ptn_ice, palb_ice, pqns_ice, pqsr_ice,   &
354         &                          pdqn_ice, pqla_ice, pdql_ice, k_limflx )
355      !!---------------------------------------------------------------------
356      !!                  ***  ROUTINE sbc_ice_lim  ***
357      !!                   
358      !! ** Purpose :   update the ice surface boundary condition by averaging and / or
359      !!                redistributing fluxes on ice categories                   
360      !!
361      !! ** Method  :   average then redistribute
362      !!
363      !! ** Action  :   
364      !!---------------------------------------------------------------------
365      INTEGER                   , INTENT(in   ) ::   k_limflx   ! =-1 do nothing; =0 average ;
366                                                                ! =1 average and redistribute ; =2 redistribute
367      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   ptn_ice    ! ice surface temperature
368      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   palb_ice   ! ice albedo
369      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pqns_ice   ! non solar flux
370      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pqsr_ice   ! net solar flux
371      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pdqn_ice   ! non solar flux sensitivity
372      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pqla_ice   ! latent heat flux
373      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pdql_ice   ! latent heat flux sensitivity
374      !
375      INTEGER  ::   jl      ! dummy loop index
376      !
377      REAL(wp), POINTER, DIMENSION(:,:) :: zalb_m    ! Mean albedo over all categories
378      REAL(wp), POINTER, DIMENSION(:,:) :: ztem_m    ! Mean temperature over all categories
379      !
380      REAL(wp), POINTER, DIMENSION(:,:) :: z_qsr_m   ! Mean solar heat flux over all categories
381      REAL(wp), POINTER, DIMENSION(:,:) :: z_qns_m   ! Mean non solar heat flux over all categories
382      REAL(wp), POINTER, DIMENSION(:,:) :: z_qla_m   ! Mean latent heat flux over all categories
383      REAL(wp), POINTER, DIMENSION(:,:) :: z_dqn_m   ! Mean d(qns)/dT over all categories
384      REAL(wp), POINTER, DIMENSION(:,:) :: z_dql_m   ! Mean d(qla)/dT over all categories
385      !!----------------------------------------------------------------------
386
387      IF( nn_timing == 1 )  CALL timing_start('ice_lim_flx')
388      !
389      !
390      SELECT CASE( k_limflx )                              !==  averaged on all ice categories  ==!
391      CASE( 0 , 1 )
392         CALL wrk_alloc( jpi,jpj, z_qsr_m, z_qns_m, z_qla_m, z_dqn_m, z_dql_m)
393         !
394         z_qns_m(:,:) = fice_ice_ave ( pqns_ice (:,:,:) )
395         z_qsr_m(:,:) = fice_ice_ave ( pqsr_ice (:,:,:) )
396         z_dqn_m(:,:) = fice_ice_ave ( pdqn_ice (:,:,:) )
397         z_qla_m(:,:) = fice_ice_ave ( pqla_ice (:,:,:) )
398         z_dql_m(:,:) = fice_ice_ave ( pdql_ice (:,:,:) )
399         DO jl = 1, jpl
400            pdqn_ice(:,:,jl) = z_dqn_m(:,:)
401            pdql_ice(:,:,jl) = z_dql_m(:,:)
402         END DO
403         !
404         DO jl = 1, jpl
405            pqns_ice(:,:,jl) = z_qns_m(:,:)
406            pqsr_ice(:,:,jl) = z_qsr_m(:,:)
407            pqla_ice(:,:,jl) = z_qla_m(:,:)
408         END DO
409         !
410         CALL wrk_dealloc( jpi,jpj, z_qsr_m, z_qns_m, z_qla_m, z_dqn_m, z_dql_m)
411      END SELECT
412
413      SELECT CASE( k_limflx )                              !==  redistribution on all ice categories  ==!
414      CASE( 1 , 2 )
415         CALL wrk_alloc( jpi,jpj, zalb_m, ztem_m )
416         !
417         zalb_m(:,:) = fice_ice_ave ( palb_ice (:,:,:) ) 
418         ztem_m(:,:) = fice_ice_ave ( ptn_ice  (:,:,:) ) 
419         DO jl = 1, jpl
420            pqns_ice(:,:,jl) = pqns_ice(:,:,jl) + pdqn_ice(:,:,jl) * (ptn_ice(:,:,jl) - ztem_m(:,:))
421            pqla_ice(:,:,jl) = pqla_ice(:,:,jl) + pdql_ice(:,:,jl) * (ptn_ice(:,:,jl) - ztem_m(:,:))
422            pqsr_ice(:,:,jl) = pqsr_ice(:,:,jl) * ( 1._wp - palb_ice(:,:,jl) ) / ( 1._wp - zalb_m(:,:) ) 
423         END DO
424         !
425         CALL wrk_dealloc( jpi,jpj, zalb_m, ztem_m )
426      END SELECT
427      !
428      IF( nn_timing == 1 )  CALL timing_stop('ice_lim_flx')
429      !
430   END SUBROUTINE ice_lim_flx
431   
432   
433   SUBROUTINE lim_ctl( kt )
434      !!-----------------------------------------------------------------------
435      !!                   ***  ROUTINE lim_ctl ***
436      !!                 
437      !! ** Purpose :   Alerts in case of model crash
438      !!-------------------------------------------------------------------
439      INTEGER, INTENT(in) ::   kt      ! ocean time step
440      INTEGER  ::   ji, jj, jk,  jl   ! dummy loop indices
441      INTEGER  ::   inb_altests       ! number of alert tests (max 20)
442      INTEGER  ::   ialert_id         ! number of the current alert
443      REAL(wp) ::   ztmelts           ! ice layer melting point
444      CHARACTER (len=30), DIMENSION(20)      ::   cl_alname   ! name of alert
445      INTEGER           , DIMENSION(20)      ::   inb_alp     ! number of alerts positive
446      !!-------------------------------------------------------------------
447
448      inb_altests = 10
449      inb_alp(:)  =  0
450
451      ! Alert if incompatible volume and concentration
452      ialert_id = 2 ! reference number of this alert
453      cl_alname(ialert_id) = ' Incompat vol and con         '    ! name of the alert
454
455      DO jl = 1, jpl
456         DO jj = 1, jpj
457            DO ji = 1, jpi
458               IF(  v_i(ji,jj,jl) /= 0._wp   .AND.   a_i(ji,jj,jl) == 0._wp   ) THEN
459                  !WRITE(numout,*) ' ALERTE 2 :   Incompatible volume and concentration '
460                  !WRITE(numout,*) ' at_i     ', at_i(ji,jj)
461                  !WRITE(numout,*) ' Point - category', ji, jj, jl
462                  !WRITE(numout,*) ' a_i *** a_i_b   ', a_i      (ji,jj,jl), a_i_b  (ji,jj,jl)
463                  !WRITE(numout,*) ' v_i *** v_i_b   ', v_i      (ji,jj,jl), v_i_b  (ji,jj,jl)
464                  !WRITE(numout,*) ' d_a_i_thd/trp   ', d_a_i_thd(ji,jj,jl), d_a_i_trp(ji,jj,jl)
465                  !WRITE(numout,*) ' d_v_i_thd/trp   ', d_v_i_thd(ji,jj,jl), d_v_i_trp(ji,jj,jl)
466                  inb_alp(ialert_id) = inb_alp(ialert_id) + 1
467               ENDIF
468            END DO
469         END DO
470      END DO
471
472      ! Alerte if very thick ice
473      ialert_id = 3 ! reference number of this alert
474      cl_alname(ialert_id) = ' Very thick ice               ' ! name of the alert
475      jl = jpl 
476      DO jj = 1, jpj
477         DO ji = 1, jpi
478            IF(   ht_i(ji,jj,jl)  >  50._wp   ) THEN
479               !CALL lim_prt_state( kt, ji, jj, 2, ' ALERTE 3 :   Very thick ice ' )
480               inb_alp(ialert_id) = inb_alp(ialert_id) + 1
481            ENDIF
482         END DO
483      END DO
484
485      ! Alert if very fast ice
486      ialert_id = 4 ! reference number of this alert
487      cl_alname(ialert_id) = ' Very fast ice               ' ! name of the alert
488      DO jj = 1, jpj
489         DO ji = 1, jpi
490            IF(   MAX( ABS( u_ice(ji,jj) ), ABS( v_ice(ji,jj) ) ) > 1.5  .AND.  &
491               &  at_i(ji,jj) > 0._wp   ) THEN
492               !CALL lim_prt_state( kt, ji, jj, 1, ' ALERTE 4 :   Very fast ice ' )
493               !WRITE(numout,*) ' ice strength             : ', strength(ji,jj)
494               !WRITE(numout,*) ' oceanic stress utau      : ', utau(ji,jj)
495               !WRITE(numout,*) ' oceanic stress vtau      : ', vtau(ji,jj)
496               !WRITE(numout,*) ' sea-ice stress utau_ice  : ', utau_ice(ji,jj)
497               !WRITE(numout,*) ' sea-ice stress vtau_ice  : ', vtau_ice(ji,jj)
498               !WRITE(numout,*) ' oceanic speed u          : ', u_oce(ji,jj)
499               !WRITE(numout,*) ' oceanic speed v          : ', v_oce(ji,jj)
500               !WRITE(numout,*) ' sst                      : ', sst_m(ji,jj)
501               !WRITE(numout,*) ' sss                      : ', sss_m(ji,jj)
502               !WRITE(numout,*)
503               inb_alp(ialert_id) = inb_alp(ialert_id) + 1
504            ENDIF
505         END DO
506      END DO
507
508      ! Alert if there is ice on continents
509      ialert_id = 6 ! reference number of this alert
510      cl_alname(ialert_id) = ' Ice on continents           ' ! name of the alert
511      DO jj = 1, jpj
512         DO ji = 1, jpi
513            IF(   tms(ji,jj) <= 0._wp   .AND.   at_i(ji,jj) > 0._wp   ) THEN 
514               !CALL lim_prt_state( kt, ji, jj, 1, ' ALERTE 6 :   Ice on continents ' )
515               !WRITE(numout,*) ' masks s, u, v        : ', tms(ji,jj), tmu(ji,jj), tmv(ji,jj)
516               !WRITE(numout,*) ' sst                  : ', sst_m(ji,jj)
517               !WRITE(numout,*) ' sss                  : ', sss_m(ji,jj)
518               !WRITE(numout,*) ' at_i(ji,jj)          : ', at_i(ji,jj)
519               !WRITE(numout,*) ' v_ice(ji,jj)         : ', v_ice(ji,jj)
520               !WRITE(numout,*) ' v_ice(ji,jj-1)       : ', v_ice(ji,jj-1)
521               !WRITE(numout,*) ' u_ice(ji-1,jj)       : ', u_ice(ji-1,jj)
522               !WRITE(numout,*) ' u_ice(ji,jj)         : ', v_ice(ji,jj)
523               !
524               inb_alp(ialert_id) = inb_alp(ialert_id) + 1
525            ENDIF
526         END DO
527      END DO
528
529!
530!     ! Alert if very fresh ice
531      ialert_id = 7 ! reference number of this alert
532      cl_alname(ialert_id) = ' Very fresh ice               ' ! name of the alert
533      DO jl = 1, jpl
534         DO jj = 1, jpj
535            DO ji = 1, jpi
536               IF( sm_i(ji,jj,jl) < 0.1 .AND. a_i(ji,jj,jl) > 0._wp ) THEN
537!                 CALL lim_prt_state(kt,ji,jj,1, ' ALERTE 7 :   Very fresh ice ' )
538!                 WRITE(numout,*) ' sst                  : ', sst_m(ji,jj)
539!                 WRITE(numout,*) ' sss                  : ', sss_m(ji,jj)
540!                 WRITE(numout,*)
541                  inb_alp(ialert_id) = inb_alp(ialert_id) + 1
542               ENDIF
543            END DO
544         END DO
545      END DO
546!
547
548!     ! Alert if too old ice
549      ialert_id = 9 ! reference number of this alert
550      cl_alname(ialert_id) = ' Very old   ice               ' ! name of the alert
551      DO jl = 1, jpl
552         DO jj = 1, jpj
553            DO ji = 1, jpi
554               IF ( ( ( ABS( o_i(ji,jj,jl) ) > rdt_ice ) .OR. &
555                      ( ABS( o_i(ji,jj,jl) ) < 0._wp) ) .AND. &
556                             ( a_i(ji,jj,jl) > 0._wp ) ) THEN
557                  !CALL lim_prt_state( kt, ji, jj, 1, ' ALERTE 9 :   Wrong ice age ')
558                  inb_alp(ialert_id) = inb_alp(ialert_id) + 1
559               ENDIF
560            END DO
561         END DO
562      END DO
563 
564      ! Alert on salt flux
565      ialert_id = 5 ! reference number of this alert
566      cl_alname(ialert_id) = ' High salt flux               ' ! name of the alert
567      DO jj = 1, jpj
568         DO ji = 1, jpi
569            IF( ABS( sfx (ji,jj) ) .GT. 1.0e-2 ) THEN  ! = 1 psu/day for 1m ocean depth
570               !CALL lim_prt_state( kt, ji, jj, 3, ' ALERTE 5 :   High salt flux ' )
571               !DO jl = 1, jpl
572                  !WRITE(numout,*) ' Category no: ', jl
573                  !WRITE(numout,*) ' a_i        : ', a_i      (ji,jj,jl) , ' a_i_b      : ', a_i_b  (ji,jj,jl)   
574                  !WRITE(numout,*) ' d_a_i_trp  : ', d_a_i_trp(ji,jj,jl) , ' d_a_i_thd  : ', d_a_i_thd(ji,jj,jl)
575                  !WRITE(numout,*) ' v_i        : ', v_i      (ji,jj,jl) , ' v_i_b      : ', v_i_b  (ji,jj,jl)   
576                  !WRITE(numout,*) ' d_v_i_trp  : ', d_v_i_trp(ji,jj,jl) , ' d_v_i_thd  : ', d_v_i_thd(ji,jj,jl)
577                  !WRITE(numout,*) ' '
578               !END DO
579               inb_alp(ialert_id) = inb_alp(ialert_id) + 1
580            ENDIF
581         END DO
582      END DO
583
584      ! Alert if qns very big
585      ialert_id = 8 ! reference number of this alert
586      cl_alname(ialert_id) = ' fnsolar very big             ' ! name of the alert
587      DO jj = 1, jpj
588         DO ji = 1, jpi
589            IF( ABS( qns(ji,jj) ) > 1500._wp .AND. at_i(ji,jj) > 0._wp ) THEN
590               !
591               !WRITE(numout,*) ' ALERTE 8 :   Very high non-solar heat flux'
592               !WRITE(numout,*) ' ji, jj    : ', ji, jj
593               !WRITE(numout,*) ' qns       : ', qns(ji,jj)
594               !WRITE(numout,*) ' sst       : ', sst_m(ji,jj)
595               !WRITE(numout,*) ' sss       : ', sss_m(ji,jj)
596               !
597               !CALL lim_prt_state( kt, ji, jj, 2, '   ')
598               inb_alp(ialert_id) = inb_alp(ialert_id) + 1
599               !
600            ENDIF
601         END DO
602      END DO
603      !+++++
604 
605      ! Alert if very warm ice
606      ialert_id = 10 ! reference number of this alert
607      cl_alname(ialert_id) = ' Very warm ice                ' ! name of the alert
608      inb_alp(ialert_id) = 0
609      DO jl = 1, jpl
610         DO jk = 1, nlay_i
611            DO jj = 1, jpj
612               DO ji = 1, jpi
613                  ztmelts    =  -tmut * s_i(ji,jj,jk,jl) + rtt
614                  IF( t_i(ji,jj,jk,jl) >= ztmelts  .AND.  v_i(ji,jj,jl) > 1.e-10   &
615                     &                             .AND.  a_i(ji,jj,jl) > 0._wp   ) THEN
616                     !WRITE(numout,*) ' ALERTE 10 :   Very warm ice'
617                     !WRITE(numout,*) ' ji, jj, jk, jl : ', ji, jj, jk, jl
618                     !WRITE(numout,*) ' t_i : ', t_i(ji,jj,jk,jl)
619                     !WRITE(numout,*) ' e_i : ', e_i(ji,jj,jk,jl)
620                     !WRITE(numout,*) ' s_i : ', s_i(ji,jj,jk,jl)
621                     !WRITE(numout,*) ' ztmelts : ', ztmelts
622                     inb_alp(ialert_id) = inb_alp(ialert_id) + 1
623                  ENDIF
624               END DO
625            END DO
626         END DO
627      END DO
628
629      ! sum of the alerts on all processors
630      IF( lk_mpp ) THEN
631         DO ialert_id = 1, inb_altests
632            CALL mpp_sum(inb_alp(ialert_id))
633         END DO
634      ENDIF
635
636      ! print alerts
637      IF( lwp ) THEN
638         ialert_id = 1                                 ! reference number of this alert
639         cl_alname(ialert_id) = ' NO alerte 1      '   ! name of the alert
640         WRITE(numout,*) ' time step ',kt
641         WRITE(numout,*) ' All alerts at the end of ice model '
642         DO ialert_id = 1, inb_altests
643            WRITE(numout,*) ialert_id, cl_alname(ialert_id)//' : ', inb_alp(ialert_id), ' times ! '
644         END DO
645      ENDIF
646     !
647   END SUBROUTINE lim_ctl
648 
649   
650   SUBROUTINE lim_prt_state( kt, ki, kj, kn, cd1 )
651      !!-----------------------------------------------------------------------
652      !!                   ***  ROUTINE lim_prt_state ***
653      !!                 
654      !! ** Purpose :   Writes global ice state on the (i,j) point
655      !!                in ocean.ouput
656      !!                3 possibilities exist
657      !!                n = 1/-1 -> simple ice state (plus Mechanical Check if -1)
658      !!                n = 2    -> exhaustive state
659      !!                n = 3    -> ice/ocean salt fluxes
660      !!
661      !! ** input   :   point coordinates (i,j)
662      !!                n : number of the option
663      !!-------------------------------------------------------------------
664      INTEGER         , INTENT(in) ::   kt            ! ocean time step
665      INTEGER         , INTENT(in) ::   ki, kj, kn    ! ocean gridpoint indices
666      CHARACTER(len=*), INTENT(in) ::   cd1           !
667      !!
668      INTEGER :: jl, ji, jj
669      !!-------------------------------------------------------------------
670
671      DO ji = mi0(ki), mi1(ki)
672         DO jj = mj0(kj), mj1(kj)
673
674            WRITE(numout,*) ' time step ',kt,' ',cd1             ! print title
675
676            !----------------
677            !  Simple state
678            !----------------
679           
680            IF ( kn == 1 .OR. kn == -1 ) THEN
681               WRITE(numout,*) ' lim_prt_state - Point : ',ji,jj
682               WRITE(numout,*) ' ~~~~~~~~~~~~~~ '
683               WRITE(numout,*) ' Simple state '
684               WRITE(numout,*) ' masks s,u,v   : ', tms(ji,jj), tmu(ji,jj), tmv(ji,jj)
685               WRITE(numout,*) ' lat - long    : ', gphit(ji,jj), glamt(ji,jj)
686               WRITE(numout,*) ' Time step     : ', numit
687               WRITE(numout,*) ' - Ice drift   '
688               WRITE(numout,*) '   ~~~~~~~~~~~ '
689               WRITE(numout,*) ' u_ice(i-1,j)  : ', u_ice(ji-1,jj)
690               WRITE(numout,*) ' u_ice(i  ,j)  : ', u_ice(ji,jj)
691               WRITE(numout,*) ' v_ice(i  ,j-1): ', v_ice(ji,jj-1)
692               WRITE(numout,*) ' v_ice(i  ,j)  : ', v_ice(ji,jj)
693               WRITE(numout,*) ' strength      : ', strength(ji,jj)
694               WRITE(numout,*)
695               WRITE(numout,*) ' - Cell values '
696               WRITE(numout,*) '   ~~~~~~~~~~~ '
697               WRITE(numout,*) ' cell area     : ', area(ji,jj)
698               WRITE(numout,*) ' at_i          : ', at_i(ji,jj)       
699               WRITE(numout,*) ' vt_i          : ', vt_i(ji,jj)       
700               WRITE(numout,*) ' vt_s          : ', vt_s(ji,jj)       
701               DO jl = 1, jpl
702                  WRITE(numout,*) ' - Category (', jl,')'
703                  WRITE(numout,*) ' a_i           : ', a_i(ji,jj,jl)
704                  WRITE(numout,*) ' ht_i          : ', ht_i(ji,jj,jl)
705                  WRITE(numout,*) ' ht_s          : ', ht_s(ji,jj,jl)
706                  WRITE(numout,*) ' v_i           : ', v_i(ji,jj,jl)
707                  WRITE(numout,*) ' v_s           : ', v_s(ji,jj,jl)
708                  WRITE(numout,*) ' e_s           : ', e_s(ji,jj,1,jl)/1.0e9
709                  WRITE(numout,*) ' e_i           : ', e_i(ji,jj,1:nlay_i,jl)/1.0e9
710                  WRITE(numout,*) ' t_su          : ', t_su(ji,jj,jl)
711                  WRITE(numout,*) ' t_snow        : ', t_s(ji,jj,1,jl)
712                  WRITE(numout,*) ' t_i           : ', t_i(ji,jj,1:nlay_i,jl)
713                  WRITE(numout,*) ' sm_i          : ', sm_i(ji,jj,jl)
714                  WRITE(numout,*) ' smv_i         : ', smv_i(ji,jj,jl)
715                  WRITE(numout,*)
716               END DO
717            ENDIF
718            IF( kn == -1 ) THEN
719               WRITE(numout,*) ' Mechanical Check ************** '
720               WRITE(numout,*) ' Check what means ice divergence '
721               WRITE(numout,*) ' Total ice concentration ', at_i (ji,jj)
722               WRITE(numout,*) ' Total lead fraction     ', ato_i(ji,jj)
723               WRITE(numout,*) ' Sum of both             ', ato_i(ji,jj) + at_i(ji,jj)
724               WRITE(numout,*) ' Sum of both minus 1     ', ato_i(ji,jj) + at_i(ji,jj) - 1.00
725            ENDIF
726           
727
728            !--------------------
729            !  Exhaustive state
730            !--------------------
731           
732            IF ( kn .EQ. 2 ) THEN
733               WRITE(numout,*) ' lim_prt_state - Point : ',ji,jj
734               WRITE(numout,*) ' ~~~~~~~~~~~~~~ '
735               WRITE(numout,*) ' Exhaustive state '
736               WRITE(numout,*) ' lat - long ', gphit(ji,jj), glamt(ji,jj)
737               WRITE(numout,*) ' Time step ', numit
738               WRITE(numout,*) 
739               WRITE(numout,*) ' - Cell values '
740               WRITE(numout,*) '   ~~~~~~~~~~~ '
741               WRITE(numout,*) ' cell area     : ', area(ji,jj)
742               WRITE(numout,*) ' at_i          : ', at_i(ji,jj)       
743               WRITE(numout,*) ' vt_i          : ', vt_i(ji,jj)       
744               WRITE(numout,*) ' vt_s          : ', vt_s(ji,jj)       
745               WRITE(numout,*) ' u_ice(i-1,j)  : ', u_ice(ji-1,jj)
746               WRITE(numout,*) ' u_ice(i  ,j)  : ', u_ice(ji,jj)
747               WRITE(numout,*) ' v_ice(i  ,j-1): ', v_ice(ji,jj-1)
748               WRITE(numout,*) ' v_ice(i  ,j)  : ', v_ice(ji,jj)
749               WRITE(numout,*) ' strength      : ', strength(ji,jj)
750               WRITE(numout,*) ' d_u_ice_dyn   : ', d_u_ice_dyn(ji,jj), ' d_v_ice_dyn   : ', d_v_ice_dyn(ji,jj)
751               WRITE(numout,*) ' u_ice_b       : ', u_ice_b(ji,jj)    , ' v_ice_b       : ', v_ice_b(ji,jj) 
752               WRITE(numout,*)
753               
754               DO jl = 1, jpl
755                  WRITE(numout,*) ' - Category (',jl,')'
756                  WRITE(numout,*) '   ~~~~~~~~         ' 
757                  WRITE(numout,*) ' ht_i       : ', ht_i(ji,jj,jl)             , ' ht_s       : ', ht_s(ji,jj,jl)
758                  WRITE(numout,*) ' t_i        : ', t_i(ji,jj,1:nlay_i,jl)
759                  WRITE(numout,*) ' t_su       : ', t_su(ji,jj,jl)             , ' t_s        : ', t_s(ji,jj,1,jl)
760                  WRITE(numout,*) ' sm_i       : ', sm_i(ji,jj,jl)             , ' o_i        : ', o_i(ji,jj,jl)
761                  WRITE(numout,*) ' a_i        : ', a_i(ji,jj,jl)              , ' a_i_b      : ', a_i_b(ji,jj,jl)   
762                  WRITE(numout,*) ' d_a_i_trp  : ', d_a_i_trp(ji,jj,jl)        , ' d_a_i_thd  : ', d_a_i_thd(ji,jj,jl) 
763                  WRITE(numout,*) ' v_i        : ', v_i(ji,jj,jl)              , ' v_i_b      : ', v_i_b(ji,jj,jl)   
764                  WRITE(numout,*) ' d_v_i_trp  : ', d_v_i_trp(ji,jj,jl)        , ' d_v_i_thd  : ', d_v_i_thd(ji,jj,jl) 
765                  WRITE(numout,*) ' v_s        : ', v_s(ji,jj,jl)              , ' v_s_b      : ', v_s_b(ji,jj,jl) 
766                  WRITE(numout,*) ' d_v_s_trp  : ', d_v_s_trp(ji,jj,jl)        , ' d_v_s_thd  : ', d_v_s_thd(ji,jj,jl)
767                  WRITE(numout,*) ' e_i1       : ', e_i(ji,jj,1,jl)/1.0e9      , ' ei1        : ', e_i_b(ji,jj,1,jl)/1.0e9 
768                  WRITE(numout,*) ' de_i1_trp  : ', d_e_i_trp(ji,jj,1,jl)/1.0e9, ' de_i1_thd  : ', d_e_i_thd(ji,jj,1,jl)/1.0e9
769                  WRITE(numout,*) ' e_i2       : ', e_i(ji,jj,2,jl)/1.0e9      , ' ei2_b      : ', e_i_b(ji,jj,2,jl)/1.0e9 
770                  WRITE(numout,*) ' de_i2_trp  : ', d_e_i_trp(ji,jj,2,jl)/1.0e9, ' de_i2_thd  : ', d_e_i_thd(ji,jj,2,jl)/1.0e9
771                  WRITE(numout,*) ' e_snow     : ', e_s(ji,jj,1,jl)            , ' e_snow_b   : ', e_s_b(ji,jj,1,jl) 
772                  WRITE(numout,*) ' d_e_s_trp  : ', d_e_s_trp(ji,jj,1,jl)      , ' d_e_s_thd  : ', d_e_s_thd(ji,jj,1,jl)
773                  WRITE(numout,*) ' smv_i      : ', smv_i(ji,jj,jl)            , ' smv_i_b    : ', smv_i_b(ji,jj,jl)   
774                  WRITE(numout,*) ' d_smv_i_trp: ', d_smv_i_trp(ji,jj,jl)      , ' d_smv_i_thd: ', d_smv_i_thd(ji,jj,jl) 
775                  WRITE(numout,*) ' oa_i       : ', oa_i(ji,jj,jl)             , ' oa_i_b     : ', oa_i_b(ji,jj,jl)
776                  WRITE(numout,*) ' d_oa_i_trp : ', d_oa_i_trp(ji,jj,jl)       , ' d_oa_i_thd : ', d_oa_i_thd(ji,jj,jl)
777               END DO !jl
778               
779               WRITE(numout,*)
780               WRITE(numout,*) ' - Heat / FW fluxes '
781               WRITE(numout,*) '   ~~~~~~~~~~~~~~~~ '
782               WRITE(numout,*) ' - Heat fluxes in and out the ice ***'
783               WRITE(numout,*) ' qsr_ini       : ', pfrld(ji,jj) * qsr(ji,jj) + SUM( a_i_b(ji,jj,:) * qsr_ice(ji,jj,:) )
784               WRITE(numout,*) ' qns_ini       : ', pfrld(ji,jj) * qns(ji,jj) + SUM( a_i_b(ji,jj,:) * qns_ice(ji,jj,:) )
785               WRITE(numout,*)
786               WRITE(numout,*) 
787               WRITE(numout,*) ' sst        : ', sst_m(ji,jj) 
788               WRITE(numout,*) ' sss        : ', sss_m(ji,jj) 
789               WRITE(numout,*) 
790               WRITE(numout,*) ' - Stresses '
791               WRITE(numout,*) '   ~~~~~~~~ '
792               WRITE(numout,*) ' utau_ice   : ', utau_ice(ji,jj) 
793               WRITE(numout,*) ' vtau_ice   : ', vtau_ice(ji,jj)
794               WRITE(numout,*) ' utau       : ', utau    (ji,jj) 
795               WRITE(numout,*) ' vtau       : ', vtau    (ji,jj)
796               WRITE(numout,*) ' oc. vel. u : ', u_oce   (ji,jj)
797               WRITE(numout,*) ' oc. vel. v : ', v_oce   (ji,jj)
798            ENDIF
799           
800            !---------------------
801            ! Salt / heat fluxes
802            !---------------------
803           
804            IF ( kn .EQ. 3 ) THEN
805               WRITE(numout,*) ' lim_prt_state - Point : ',ji,jj
806               WRITE(numout,*) ' ~~~~~~~~~~~~~~ '
807               WRITE(numout,*) ' - Salt / Heat Fluxes '
808               WRITE(numout,*) '   ~~~~~~~~~~~~~~~~ '
809               WRITE(numout,*) ' lat - long ', gphit(ji,jj), glamt(ji,jj)
810               WRITE(numout,*) ' Time step ', numit
811               WRITE(numout,*)
812               WRITE(numout,*) ' - Heat fluxes at bottom interface ***'
813               WRITE(numout,*) ' qsr       : ', qsr(ji,jj)
814               WRITE(numout,*) ' qns       : ', qns(ji,jj)
815               WRITE(numout,*)
816               WRITE(numout,*) ' hfx_mass     : ', hfx_thd(ji,jj) + hfx_dyn(ji,jj) + hfx_snw(ji,jj) + hfx_res(ji,jj)
817               WRITE(numout,*) ' hfx_in       : ', hfx_in(ji,jj)
818               WRITE(numout,*) ' hfx_out      : ', hfx_out(ji,jj)
819               WRITE(numout,*) ' dhc          : ', diag_heat_dhc(ji,jj)             
820               WRITE(numout,*)
821               WRITE(numout,*) ' hfx_dyn      : ', hfx_dyn(ji,jj)
822               WRITE(numout,*) ' hfx_thd      : ', hfx_thd(ji,jj)
823               WRITE(numout,*) ' hfx_res      : ', hfx_res(ji,jj)
824               WRITE(numout,*) ' fhtur        : ', fhtur(ji,jj) 
825               WRITE(numout,*) ' qlead        : ', qlead(ji,jj) * r1_rdtice
826               WRITE(numout,*)
827               WRITE(numout,*) ' - Salt fluxes at bottom interface ***'
828               WRITE(numout,*) ' emp       : ', emp    (ji,jj)
829               WRITE(numout,*) ' sfx       : ', sfx    (ji,jj)
830               WRITE(numout,*) ' sfx_res   : ', sfx_res(ji,jj)
831               WRITE(numout,*) ' sfx_bri   : ', sfx_bri(ji,jj)
832               WRITE(numout,*) ' sfx_dyn   : ', sfx_dyn(ji,jj)
833               WRITE(numout,*)
834               WRITE(numout,*) ' - Momentum fluxes '
835               WRITE(numout,*) ' utau      : ', utau(ji,jj) 
836               WRITE(numout,*) ' vtau      : ', vtau(ji,jj)
837            ENDIF
838            WRITE(numout,*) ' '
839            !
840         END DO
841      END DO
842      !
843   END SUBROUTINE lim_prt_state
844   
845     
846   FUNCTION fice_cell_ave ( ptab )
847      !!--------------------------------------------------------------------------
848      !! * Compute average over categories, for grid cell (ice covered and free ocean)
849      !!--------------------------------------------------------------------------
850      REAL (wp), DIMENSION (jpi,jpj) :: fice_cell_ave
851      REAL (wp), DIMENSION (jpi,jpj,jpl), INTENT (in) :: ptab
852      INTEGER :: jl ! Dummy loop index
853     
854      fice_cell_ave (:,:) = 0.0_wp
855     
856      DO jl = 1, jpl
857         fice_cell_ave (:,:) = fice_cell_ave (:,:) &
858            &                  + a_i (:,:,jl) * ptab (:,:,jl)
859      END DO
860     
861   END FUNCTION fice_cell_ave
862   
863   
864   FUNCTION fice_ice_ave ( ptab )
865      !!--------------------------------------------------------------------------
866      !! * Compute average over categories, for ice covered part of grid cell
867      !!--------------------------------------------------------------------------
868      REAL (kind=wp), DIMENSION (jpi,jpj) :: fice_ice_ave
869      REAL (kind=wp), DIMENSION (jpi,jpj,jpl), INTENT(in) :: ptab
870
871      fice_ice_ave (:,:) = 0.0_wp
872      WHERE ( at_i (:,:) .GT. 0.0_wp ) fice_ice_ave (:,:) = fice_cell_ave ( ptab (:,:,:)) / at_i (:,:)
873
874   END FUNCTION fice_ice_ave
875
876
877#else
878   !!----------------------------------------------------------------------
879   !!   Default option           Dummy module      NO LIM 3.0 sea-ice model
880   !!----------------------------------------------------------------------
881CONTAINS
882   SUBROUTINE sbc_ice_lim ( kt, kblk )     ! Dummy routine
883      WRITE(*,*) 'sbc_ice_lim: You should not have seen this print! error?', kt, kblk
884   END SUBROUTINE sbc_ice_lim
885#endif
886
887   !!======================================================================
888END MODULE sbcice_lim
Note: See TracBrowser for help on using the repository browser.