source: CONFIG/UNIFORM/v6/IPSLCM6/SOURCES/NEMO/sbcice_lim.F90 @ 2131

Last change on this file since 2131 was 2113, checked in by omamce, 11 years ago

O.M.

  • Add coupled interface for LIM3
  • Closea : spread Black Sea outflow on several points
  • stpctl : prevent Salinity to be below 0.1 PSS
File size: 40.8 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   !!----------------------------------------------------------------------
14#if defined key_lim3
15   !!----------------------------------------------------------------------
16   !!   'key_lim3' :                                  LIM 3.0 sea-ice model
17   !!----------------------------------------------------------------------
18   !!   sbc_ice_lim  : sea-ice model time-stepping and update ocean sbc over ice-covered area
19   !!   lim_ctl       : alerts in case of ice model crash
20   !!   lim_prt_state : ice control print at a given grid point
21   !!----------------------------------------------------------------------
22   USE oce             ! ocean dynamics and tracers
23   USE dom_oce         ! ocean space and time domain
24   USE par_ice         ! sea-ice parameters
25   USE ice             ! LIM-3: ice variables
26   USE iceini          ! LIM-3: ice initialisation
27   USE dom_ice         ! LIM-3: ice domain
28
29   USE sbc_oce         ! Surface boundary condition: ocean fields
30   USE sbc_ice         ! Surface boundary condition: ice   fields
31   USE sbcblk_core     ! Surface boundary condition: CORE bulk
32   USE sbcblk_clio     ! Surface boundary condition: CLIO bulk
33   USE sbccpl          ! Surface boundary condition: coupled interface
34   USE albedo          ! ocean & ice albedo
35
36   USE phycst          ! Define parameters for the routines
37   USE eosbn2          ! equation of state
38   USE limdyn          ! Ice dynamics
39   USE limtrp          ! Ice transport
40   USE limthd          ! Ice thermodynamics
41   USE limitd_th       ! Thermodynamics on ice thickness distribution
42   USE limitd_me       ! Mechanics on ice thickness distribution
43   USE limsbc          ! sea surface boundary condition
44   USE limdia          ! Ice diagnostics
45   USE limwri          ! Ice outputs
46   USE limrst          ! Ice restarts
47   USE limupdate       ! update of global variables
48   USE limvar          ! Ice variables switch
49
50   USE c1d             ! 1D vertical configuration
51   USE lbclnk          ! lateral boundary condition - MPP link
52   USE lib_mpp         ! MPP library
53   USE wrk_nemo        ! work arrays
54   USE iom             ! I/O manager library
55   USE in_out_manager  ! I/O manager
56   USE prtctl          ! Print control
57
58   IMPLICIT NONE
59   PRIVATE
60
61   REAL (kind=wp), SAVE :: rdt_ice_dyn, rdt_ice_thermo
62
63   PUBLIC sbc_ice_lim  ! routine called by sbcmod.F90
64   
65   !! * Substitutions
66#  include "domzgr_substitute.h90"
67#  include "vectopt_loop_substitute.h90"
68   !!----------------------------------------------------------------------
69   !! NEMO/OPA 4.0 , UCL NEMO Consortium (2011)
70   !! $Id: sbcice_lim.F90 3625 2012-11-21 13:19:18Z acc $
71   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
72   !!----------------------------------------------------------------------
73CONTAINS
74
75   FUNCTION fice_cell_ave ( ptab)
76      !!--------------------------------------------------------------------------
77      !! * Compute average over categories, for grid cell (ice covered and free ocean)
78      !!--------------------------------------------------------------------------
79      REAL (wp), DIMENSION (jpi,jpj) :: fice_cell_ave
80      REAL (wp), DIMENSION (jpi,jpj,jpl), INTENT (in) :: ptab
81      INTEGER :: jl ! Dummy loop index
82     
83      fice_cell_ave (:,:) = 0.0_wp
84     
85      DO jl = 1, jpl
86         fice_cell_ave (:,:) = fice_cell_ave (:,:) &
87            &                  + a_i (:,:,jl) * ptab (:,:,jl)
88      END DO
89     
90   END FUNCTION fice_cell_ave
91   
92   FUNCTION fice_ice_ave ( ptab)
93      !!--------------------------------------------------------------------------
94      !! * Compute average over categories, for ice covered part of grid cell
95      !!--------------------------------------------------------------------------
96      REAL (kind=wp), DIMENSION (jpi,jpj) :: fice_ice_ave
97      REAL (kind=wp), DIMENSION (jpi,jpj,jpl), INTENT(in) :: ptab
98
99      fice_ice_ave (:,:) = 0.0_wp
100      WHERE ( at_i (:,:) .GT. 0.0_wp ) fice_ice_ave (:,:) = fice_cell_ave ( ptab (:,:,:)) / at_i (:,:)
101
102   END FUNCTION fice_ice_ave
103
104   !!======================================================================
105
106   SUBROUTINE sbc_ice_lim( kt, kblk )
107      !!---------------------------------------------------------------------
108      !!                  ***  ROUTINE sbc_ice_lim  ***
109      !!                   
110      !! ** Purpose :   update the ocean surface boundary condition via the
111      !!                Louvain la Neuve Sea Ice Model time stepping
112      !!
113      !! ** Method  :   ice model time stepping
114      !!              - call the ice dynamics routine
115      !!              - call the ice advection/diffusion routine
116      !!              - call the ice thermodynamics routine
117      !!              - call the routine that computes mass and
118      !!                heat fluxes at the ice/ocean interface
119      !!              - save the outputs
120      !!              - save the outputs for restart when necessary
121      !!
122      !! ** Action  : - time evolution of the LIM sea-ice model
123      !!              - update all sbc variables below sea-ice:
124      !!                utau, vtau, taum, wndm, qns , qsr, emp , sfx
125      !!---------------------------------------------------------------------
126      INTEGER, INTENT(in) ::   kt      ! ocean time step
127      INTEGER, INTENT(in) ::   kblk    ! type of bulk (=3 CLIO, =4 CORE)
128      !!
129      INTEGER  ::   jl      ! dummy loop index
130      REAL(wp) ::   zcoef   ! local scalar
131      REAL(wp), POINTER, DIMENSION(:,:,:)   ::   zalb_ice_os, zalb_ice_cs  ! albedo of the ice under overcast/clear sky
132      REAL(wp), POINTER, DIMENSION(:,:,:)   ::   zalb_ice      ! mean albedo of ice (for coupled)
133
134      REAL(wp), POINTER, DIMENSION(:,:) :: zalb_ice_all    ! Mean albedo over all categories
135      REAL(wp), POINTER, DIMENSION(:,:) :: ztem_ice_all    ! Mean temperature over all categories
136     
137      REAL(wp), POINTER, DIMENSION(:,:) :: z_qsr_ice_all   ! Mean solar heat flux over all categories
138      REAL(wp), POINTER, DIMENSION(:,:) :: z_qns_ice_all   ! Mean non solar heat flux over all categories
139      REAL(wp), POINTER, DIMENSION(:,:) :: z_qla_ice_all   ! Mean latent heat flux over all categories
140      REAL(wp), POINTER, DIMENSION(:,:) :: z_dqns_ice_all  ! Mean d(qns)/dT over all categories
141      REAL(wp), POINTER, DIMENSION(:,:) :: z_dqla_ice_all  ! Mean d(qla)/dT over all categories
142      !!----------------------------------------------------------------------
143
144      !- O.M. : why do we allocate all these arrays even when MOD( kt-1, nn_fsbc ) /= 0 ?????
145
146      CALL wrk_alloc( jpi,jpj,jpl, zalb_ice_os, zalb_ice_cs )
147
148      IF ( ln_cpl .OR. ln_iceflx_ave .OR. ln_iceflx_linear ) THEN
149         CALL wrk_alloc( jpi,jpj,jpl, zalb_ice)
150      END IF
151      IF ( ln_iceflx_ave .OR. ln_iceflx_linear ) THEN
152         CALL wrk_alloc( jpi,jpj, ztem_ice_all, zalb_ice_all, z_qsr_ice_all, z_qns_ice_all, z_qla_ice_all, z_dqns_ice_all, z_dqla_ice_all)
153      ENDIF
154
155
156      IF( kt == nit000 ) THEN
157         IF(lwp) WRITE(numout,*)
158         IF(lwp) WRITE(numout,*) 'sbc_ice_lim : update ocean surface boudary condition' 
159         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   via Louvain la Neuve Ice Model (LIM-3) time stepping'
160         !
161         CALL ice_init
162         !
163         rdt_ice_dyn    = rdt_ice * REAL ( nn_ice_dyn, wp)
164         rdt_ice_thermo = rdt_ice
165         !
166         IF(lwp) WRITE(numout,*) 'nn_ice_dyn     : ', nn_ice_dyn
167         IF(lwp) WRITE(numout,*) 'rdt_ice        : ', rdt_ice
168         IF(lwp) WRITE(numout,*) 'rdt_ice_thermo : ', rdt_ice_thermo
169         IF(lwp) WRITE(numout,*) 'rdt_ice_dyn    : ', rdt_ice_dyn
170         !
171         IF( ln_nicep ) THEN      ! control print at a given point
172            jiindx = 44   ;   jjindx = 140
173            WRITE(numout,*) ' The debugging point is : jiindx : ',jiindx, ' jjindx : ',jjindx
174         ENDIF
175      ENDIF
176
177      !                                        !----------------------!
178      IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN     !  Ice time-step only  !
179         !                                     !----------------------!
180         !                                           !  Bulk Formulea !
181         !                                           !----------------!
182         !
183         u_oce(:,:) = ssu_m(:,:)                     ! mean surface ocean current at ice velocity point
184         v_oce(:,:) = ssv_m(:,:)                     ! (C-grid dynamics :  U- & V-points as the ocean)
185         !
186         t_bo(:,:) = tfreez( sss_m ) +  rt0          ! masked sea surface freezing temperature [Kelvin]
187         !                                           ! (set to rt0 over land)
188         CALL albedo_ice( t_su, ht_i, ht_s, zalb_ice_cs, zalb_ice_os )  ! ... ice albedo
189
190         DO jl = 1, jpl
191            t_su(:,:,jl) = t_su(:,:,jl) +  rt0 * ( 1. - tmask(:,:,1) )
192         END DO
193
194         IF ( ln_cpl ) zalb_ice (:,:,:) = 0.5 * ( zalb_ice_cs (:,:,:) +  zalb_ice_os (:,:,:) )
195         
196         IF ( ln_iceflx_ave .OR. ln_iceflx_linear ) THEN
197            !
198            ! Compute mean albedo and temperature
199            zalb_ice_all (:,:) = fice_ice_ave ( zalb_ice (:,:,:) ) 
200            ztem_ice_all (:,:) = fice_ice_ave ( tn_ice   (:,:,:) ) 
201            !
202         ENDIF
203                                                     ! Bulk formulea - provides the following fields:
204         ! utau_ice, vtau_ice : surface ice stress                     (U- & V-points)   [N/m2]
205         ! qsr_ice , qns_ice  : solar & non solar heat flux over ice   (T-point)         [W/m2]
206         ! qla_ice            : latent heat flux over ice              (T-point)         [W/m2]
207         ! dqns_ice, dqla_ice : non solar & latent heat sensistivity   (T-point)         [W/m2]
208         ! tprecip , sprecip  : total & solid precipitation            (T-point)         [Kg/m2/s]
209         ! fr1_i0  , fr2_i0   : 1sr & 2nd fraction of qsr penetration in ice             [%]
210         !
211         SELECT CASE( kblk )
212         CASE( 3 )                                       ! CLIO bulk formulation
213            CALL blk_ice_clio( t_su , zalb_ice_cs, zalb_ice_os,                           &
214               &                      utau_ice   , vtau_ice   , qns_ice   , qsr_ice   ,   &
215               &                      qla_ice    , dqns_ice   , dqla_ice  ,               &
216               &                      tprecip    , sprecip    ,                           &
217               &                      fr1_i0     , fr2_i0     , cp_ice_msh, jpl  )
218            !         
219         CASE( 4 )                                       ! CORE bulk formulation
220            CALL blk_ice_core( t_su , u_ice     , v_ice     , zalb_ice_cs,               &
221               &                      utau_ice  , vtau_ice  , qns_ice    , qsr_ice   ,   &
222               &                      qla_ice   , dqns_ice  , dqla_ice   ,               &
223               &                      tprecip   , sprecip   ,                            &
224               &                      fr1_i0    , fr2_i0    , cp_ice_msh, jpl  )
225            !
226         CASE ( 5 )
227            zalb_ice (:,:,:) = 0.5 * ( zalb_ice_cs (:,:,:) +  zalb_ice_os (:,:,:) )
228           
229            CALL sbc_cpl_ice_tau( utau_ice , vtau_ice )
230
231            CALL sbc_cpl_ice_flx( p_frld=ato_i, palbi=zalb_ice, psst=sst_m, pist=tn_ice    )
232
233            ! Latent heat flux is forced to 0 in coupled :
234            !  it is included in qns (non-solar heat flux)
235            qla_ice  (:,:,:) = 0.0e0_wp
236            dqla_ice (:,:,:) = 0.0e0_wp
237            !
238         END SELECT
239
240         ! Average over all categories
241         IF ( ln_iceflx_ave .OR. ln_iceflx_linear ) THEN
242
243            z_qns_ice_all  (:,:) = fice_ice_ave ( qns_ice  (:,:,:) )
244            z_qsr_ice_all  (:,:) = fice_ice_ave ( qsr_ice  (:,:,:) )
245            z_dqns_ice_all (:,:) = fice_ice_ave ( dqns_ice (:,:,:) )
246            z_qla_ice_all  (:,:) = fice_ice_ave ( qla_ice  (:,:,:) )
247            z_dqla_ice_all (:,:) = fice_ice_ave ( dqla_ice (:,:,:) )
248
249            DO jl = 1, jpl
250               dqns_ice (:,:,jl) = z_dqns_ice_all (:,:)
251               dqla_ice (:,:,jl) = z_dqla_ice_all (:,:)
252            END DO
253            !
254            IF ( ln_iceflx_ave ) THEN
255               DO jl = 1, jpl
256                  qns_ice  (:,:,jl) = z_qns_ice_all  (:,:)
257                  qsr_ice  (:,:,jl) = z_qsr_ice_all  (:,:)
258                  qla_ice  (:,:,jl) = z_qla_ice_all  (:,:)
259               END DO
260            END IF
261            !
262            IF ( ln_iceflx_linear ) THEN
263               DO jl = 1, jpl
264                  qns_ice  (:,:,jl) = z_qns_ice_all(:,:) + z_dqns_ice_all(:,:) * (tn_ice(:,:,jl) - ztem_ice_all(:,:))
265                  qla_ice  (:,:,jl) = z_qla_ice_all(:,:) + z_dqla_ice_all(:,:) * (tn_ice(:,:,jl) - ztem_ice_all(:,:))
266                  qsr_ice  (:,:,jl) = (1.0e0_wp-zalb_ice(:,:,jl)) / (1.0e0_wp-zalb_ice_all(:,:)) * z_qsr_ice_all(:,:)
267               END DO
268            END IF
269         END IF
270
271         !                                           !----------------------!
272         !                                           ! LIM-3  time-stepping !
273         !                                           !----------------------!
274         !
275         numit = numit + nn_fsbc                     ! Ice model time step
276         !
277         !                                           ! Store previous ice values
278!!gm : remark   old_...   should becomes ...b  as tn versus tb 
279         old_a_i  (:,:,:)   = a_i  (:,:,:)     ! ice area
280         old_e_i  (:,:,:,:) = e_i  (:,:,:,:)   ! ice thermal energy
281         old_v_i  (:,:,:)   = v_i  (:,:,:)     ! ice volume
282         old_v_s  (:,:,:)   = v_s  (:,:,:)     ! snow volume
283         old_e_s  (:,:,:,:) = e_s  (:,:,:,:)   ! snow thermal energy
284         old_smv_i(:,:,:)   = smv_i(:,:,:)     ! salt content
285         old_oa_i (:,:,:)   = oa_i (:,:,:)     ! areal age content
286
287         !                                           ! intialisation to zero    !!gm is it truly necessary ???
288         d_a_i_thd  (:,:,:)   = 0._wp   ;   d_a_i_trp  (:,:,:)   = 0._wp
289         d_v_i_thd  (:,:,:)   = 0._wp   ;   d_v_i_trp  (:,:,:)   = 0._wp
290         d_e_i_thd  (:,:,:,:) = 0._wp   ;   d_e_i_trp  (:,:,:,:) = 0._wp
291         d_v_s_thd  (:,:,:)   = 0._wp   ;   d_v_s_trp  (:,:,:)   = 0._wp
292         d_e_s_thd  (:,:,:,:) = 0._wp   ;   d_e_s_trp  (:,:,:,:) = 0._wp
293         d_smv_i_thd(:,:,:)   = 0._wp   ;   d_smv_i_trp(:,:,:)   = 0._wp
294         d_oa_i_thd (:,:,:)   = 0._wp   ;   d_oa_i_trp (:,:,:)   = 0._wp
295         !
296         sfx    (:,:) = 0._wp
297         sfx_bri(:,:) = 0._wp   ;   sfx_mec  (:,:) = 0._wp   ;   sfx_res  (:,:) = 0._wp
298         fhbri  (:,:) = 0._wp   ;   fheat_mec(:,:) = 0._wp   ;   fheat_res(:,:) = 0._wp
299         fhmec  (:,:) = 0._wp   ;   
300         fmmec  (:,:) = 0._wp     
301         focea2D(:,:) = 0._wp
302         fsup2D (:,:) = 0._wp
303         !
304         diag_sni_gr(:,:) = 0._wp   ;   diag_lat_gr(:,:) = 0._wp
305         diag_bot_gr(:,:) = 0._wp   ;   diag_dyn_gr(:,:) = 0._wp
306         diag_bot_me(:,:) = 0._wp   ;   diag_sur_me(:,:) = 0._wp
307         ! dynamical invariants
308         delta_i(:,:) = 0._wp       ;   divu_i(:,:) = 0._wp       ;   shear_i(:,:) = 0._wp
309
310                          CALL lim_rst_opn( kt )     ! Open Ice restart file
311         !
312         IF( ln_nicep )   CALL lim_prt_state( jiindx, jjindx, 1, ' - Beginning the time step - ' )   ! control print
313         !
314         IF( .NOT. lk_c1d ) THEN                     ! Ice dynamics & transport (except in 1D case)
315            ! - Test OM, call ice dynamics every nn_ice_dyn ice time steps
316            IF( MOD( kt-1, nn_fsbc*nn_ice_dyn ) == 0 ) THEN
317                          rdt_ice = rdt_ice_dyn
318                          CALL lim_dyn( kt )              ! Ice dynamics    ( rheology/dynamics )
319                          CALL lim_trp( kt )              ! Ice transport   ( Advection/diffusion )
320                          CALL lim_var_agg(1)             ! aggregate categories, requested
321                          CALL lim_var_glo2eqv            ! equivalent variables, requested for rafting
322         IF( ln_nicep )   CALL lim_prt_state( jiindx, jjindx,-1, ' - ice dyn & trp - ' )   ! control print
323                          CALL lim_itd_me                 ! Mechanical redistribution ! (ridging/rafting)
324            ENDIF
325         END IF
326                          rdt_ice = rdt_ice_thermo
327         !                                           ! Ice thermodynamics
328                          CALL lim_var_glo2eqv            ! equivalent variables
329                          CALL lim_var_agg(1)             ! aggregate ice categories
330                          CALL lim_var_bv                 ! bulk brine volume (diag)
331                          CALL lim_thd( kt )              ! Ice thermodynamics
332                          zcoef = rdt_ice /rday           !  Ice natural aging
333                          oa_i(:,:,:) = oa_i(:,:,:) + a_i(:,:,:) * zcoef
334                          CALL lim_var_glo2eqv            ! this CALL is maybe not necessary (Martin)
335         IF( ln_nicep )   CALL lim_prt_state( jiindx, jjindx, 1, ' - ice thermodyn. - ' )   ! control print
336                          CALL lim_itd_th( kt )           !  Remap ice categories, lateral accretion  !
337         !
338         !                                           ! Global variables update
339                          CALL lim_var_agg( 1 )           ! requested by limupdate
340                          CALL lim_update                 ! Global variables update
341                          CALL lim_var_glo2eqv            ! equivalent variables (outputs)
342                          CALL lim_var_agg(2)             ! aggregate ice thickness categories
343         IF( ln_nicep )   CALL lim_prt_state( jiindx, jjindx, 2, ' - Final state - ' )   ! control print
344         !
345                          CALL lim_sbc_flx( kt )     ! Update surface ocean mass, heat and salt fluxes
346         !
347         IF( ln_nicep )   CALL lim_prt_state( jiindx, jjindx, 3, ' - Final state lim_sbc - ' )   ! control print
348         !
349         !                                           ! Diagnostics and outputs
350         IF( ( MOD( kt+nn_fsbc-1, ninfo ) == 0 .OR. ntmoy == 1 ) .AND. .NOT. lk_mpp )   &
351            &             CALL lim_dia 
352                          CALL lim_wri( 1  )              ! Ice outputs
353         IF( lrst_ice )   CALL lim_rst_write( kt )        ! Ice restart file
354                          CALL lim_var_glo2eqv            ! ???
355         !
356         IF( ln_nicep )   CALL lim_ctl               ! alerts in case of model crash
357         !
358      ENDIF                                    ! End sea-ice time step only
359
360      !                                        !--------------------------!
361      !                                        !  at all ocean time step  !
362      !                                        !--------------------------!
363      !                                               
364      !                                              ! Update surface ocean stresses (only in ice-dynamic case)
365      !                                                   ! otherwise the atm.-ocean stresses are used everywhere
366      IF( ln_limdyn )     CALL lim_sbc_tau( kt, ub(:,:,1), vb(:,:,1) )  ! using before instantaneous surf. currents
367     
368!!gm   remark, the ocean-ice stress is not saved in ice diag call above .....  find a solution!!!
369      !
370      CALL wrk_dealloc( jpi,jpj,jpl, zalb_ice_os, zalb_ice_cs )
371      IF ( ln_cpl .OR. ln_iceflx_ave .OR. ln_iceflx_linear ) THEN
372         CALL wrk_dealloc( jpi,jpj,jpl, zalb_ice)
373      END IF
374      IF ( ln_iceflx_ave .OR. ln_iceflx_linear ) THEN
375         CALL wrk_dealloc( jpi,jpj, ztem_ice_all, zalb_ice_all, z_qsr_ice_all, z_qns_ice_all, z_qla_ice_all, z_dqns_ice_all, z_dqla_ice_all)
376      ENDIF
377      !
378   END SUBROUTINE sbc_ice_lim
379
380
381   SUBROUTINE lim_ctl
382      !!-----------------------------------------------------------------------
383      !!                   ***  ROUTINE lim_ctl ***
384      !!                 
385      !! ** Purpose :   Alerts in case of model crash
386      !!-------------------------------------------------------------------
387      INTEGER  ::   ji, jj, jk,  jl   ! dummy loop indices
388      INTEGER  ::   inb_altests       ! number of alert tests (max 20)
389      INTEGER  ::   ialert_id         ! number of the current alert
390      REAL(wp) ::   ztmelts           ! ice layer melting point
391      CHARACTER (len=30), DIMENSION(20)      ::   cl_alname   ! name of alert
392      INTEGER           , DIMENSION(20)      ::   inb_alp     ! number of alerts positive
393      !!-------------------------------------------------------------------
394
395      inb_altests = 10
396      inb_alp(:)  =  0
397
398      ! Alert if incompatible volume and concentration
399      ialert_id = 2 ! reference number of this alert
400      cl_alname(ialert_id) = ' Incompat vol and con         '    ! name of the alert
401
402      DO jl = 1, jpl
403         DO jj = 1, jpj
404            DO ji = 1, jpi
405               IF(  v_i(ji,jj,jl) /= 0._wp   .AND.   a_i(ji,jj,jl) == 0._wp   ) THEN
406                  WRITE(numout,*) ' ALERTE 2 :   Incompatible volume and concentration '
407                  WRITE(numout,*) ' at_i     ', at_i(ji,jj)
408                  WRITE(numout,*) ' Point - category', ji, jj, jl
409                  WRITE(numout,*) ' a_i *** a_i_old ', a_i      (ji,jj,jl), old_a_i  (ji,jj,jl)
410                  WRITE(numout,*) ' v_i *** v_i_old ', v_i      (ji,jj,jl), old_v_i  (ji,jj,jl)
411                  WRITE(numout,*) ' d_a_i_thd/trp   ', d_a_i_thd(ji,jj,jl), d_a_i_trp(ji,jj,jl)
412                  WRITE(numout,*) ' d_v_i_thd/trp   ', d_v_i_thd(ji,jj,jl), d_v_i_trp(ji,jj,jl)
413                  inb_alp(ialert_id) = inb_alp(ialert_id) + 1
414               ENDIF
415            END DO
416         END DO
417      END DO
418
419      ! Alerte if very thick ice
420      ialert_id = 3 ! reference number of this alert
421      cl_alname(ialert_id) = ' Very thick ice               ' ! name of the alert
422      jl = jpl 
423      DO jj = 1, jpj
424         DO ji = 1, jpi
425            IF(   ht_i(ji,jj,jl)  >  50._wp   ) THEN
426               CALL lim_prt_state( ji, jj, 2, ' ALERTE 3 :   Very thick ice ' )
427               inb_alp(ialert_id) = inb_alp(ialert_id) + 1
428            ENDIF
429         END DO
430      END DO
431
432      ! Alert if very fast ice
433      ialert_id = 4 ! reference number of this alert
434      cl_alname(ialert_id) = ' Very fast ice               ' ! name of the alert
435      DO jj = 1, jpj
436         DO ji = 1, jpi
437            IF(   MAX( ABS( u_ice(ji,jj) ), ABS( v_ice(ji,jj) ) ) > 0.5  .AND.  &
438               &  at_i(ji,jj) > 0._wp   ) THEN
439               CALL lim_prt_state( ji, jj, 1, ' ALERTE 4 :   Very fast ice ' )
440               WRITE(numout,*) ' ice strength             : ', strength(ji,jj)
441               WRITE(numout,*) ' oceanic stress utau      : ', utau(ji,jj) 
442               WRITE(numout,*) ' oceanic stress vtau      : ', vtau(ji,jj)
443               WRITE(numout,*) ' sea-ice stress utau_ice  : ', utau_ice(ji,jj) 
444               WRITE(numout,*) ' sea-ice stress vtau_ice  : ', vtau_ice(ji,jj)
445               WRITE(numout,*) ' oceanic speed u          : ', u_oce(ji,jj)
446               WRITE(numout,*) ' oceanic speed v          : ', v_oce(ji,jj)
447               WRITE(numout,*) ' sst                      : ', sst_m(ji,jj)
448               WRITE(numout,*) ' sss                      : ', sss_m(ji,jj)
449               WRITE(numout,*) 
450               inb_alp(ialert_id) = inb_alp(ialert_id) + 1
451            ENDIF
452         END DO
453      END DO
454
455      ! Alert if there is ice on continents
456      ialert_id = 6 ! reference number of this alert
457      cl_alname(ialert_id) = ' Ice on continents           ' ! name of the alert
458      DO jj = 1, jpj
459         DO ji = 1, jpi
460            IF(   tms(ji,jj) <= 0._wp   .AND.   at_i(ji,jj) > 0._wp   ) THEN
461               CALL lim_prt_state( ji, jj, 1, ' ALERTE 6 :   Ice on continents ' )
462               WRITE(numout,*) ' masks s, u, v        : ', tms(ji,jj), tmu(ji,jj), tmv(ji,jj) 
463               WRITE(numout,*) ' sst                  : ', sst_m(ji,jj)
464               WRITE(numout,*) ' sss                  : ', sss_m(ji,jj)
465               WRITE(numout,*) ' at_i(ji,jj)          : ', at_i(ji,jj)
466               WRITE(numout,*) ' v_ice(ji,jj)         : ', v_ice(ji,jj)
467               WRITE(numout,*) ' v_ice(ji,jj-1)       : ', v_ice(ji,jj-1)
468               WRITE(numout,*) ' u_ice(ji-1,jj)       : ', u_ice(ji-1,jj)
469               WRITE(numout,*) ' u_ice(ji,jj)         : ', v_ice(ji,jj)
470               !
471               inb_alp(ialert_id) = inb_alp(ialert_id) + 1
472            ENDIF
473         END DO
474      END DO
475
476!
477!     ! Alert if very fresh ice
478      ialert_id = 7 ! reference number of this alert
479      cl_alname(ialert_id) = ' Very fresh ice               ' ! name of the alert
480      DO jl = 1, jpl
481         DO jj = 1, jpj
482            DO ji = 1, jpi
483!!gm  test twice sm_i ...  ????  bug?
484               IF( ( ( ABS( sm_i(ji,jj,jl) ) < 0.5 )   .OR.  &
485                     ( ABS( sm_i(ji,jj,jl) ) < 0.5 ) ) .AND. &
486                             ( a_i(ji,jj,jl) > 0._wp ) ) THEN
487!                 CALL lim_prt_state(ji,jj,1, ' ALERTE 7 :   Very fresh ice ' )
488!                 WRITE(numout,*) ' sst                  : ', sst_m(ji,jj)
489!                 WRITE(numout,*) ' sss                  : ', sss_m(ji,jj)
490!                 WRITE(numout,*) ' s_i_newice           : ', s_i_newice(ji,jj,1:jpl)
491!                 WRITE(numout,*)
492                  inb_alp(ialert_id) = inb_alp(ialert_id) + 1
493               ENDIF
494            END DO
495         END DO
496      END DO
497!
498
499!     ! Alert if too old ice
500      ialert_id = 9 ! reference number of this alert
501      cl_alname(ialert_id) = ' Very old   ice               ' ! name of the alert
502      DO jl = 1, jpl
503         DO jj = 1, jpj
504            DO ji = 1, jpi
505               IF ( ( ( ABS( o_i(ji,jj,jl) ) > rdt_ice ) .OR. &
506                      ( ABS( o_i(ji,jj,jl) ) < 0._wp) ) .AND. &
507                             ( a_i(ji,jj,jl) > 0._wp ) ) THEN
508                  CALL lim_prt_state( ji, jj, 1, ' ALERTE 9 :   Wrong ice age ')
509                  inb_alp(ialert_id) = inb_alp(ialert_id) + 1
510               ENDIF
511            END DO
512         END DO
513      END DO
514 
515      ! Alert on salt flux
516      ialert_id = 5 ! reference number of this alert
517      cl_alname(ialert_id) = ' High salt flux               ' ! name of the alert
518      DO jj = 1, jpj
519         DO ji = 1, jpi
520            IF( ABS( sfx (ji,jj) ) .GT. 1.0e-2 ) THEN
521               CALL lim_prt_state( ji, jj, 3, ' ALERTE 5 :   High salt flux ' )
522               DO jl = 1, jpl
523                  WRITE(numout,*) ' Category no: ', jl
524                  WRITE(numout,*) ' a_i        : ', a_i      (ji,jj,jl) , ' old_a_i    : ', old_a_i  (ji,jj,jl)   
525                  WRITE(numout,*) ' d_a_i_trp  : ', d_a_i_trp(ji,jj,jl) , ' d_a_i_thd  : ', d_a_i_thd(ji,jj,jl) 
526                  WRITE(numout,*) ' v_i        : ', v_i      (ji,jj,jl) , ' old_v_i    : ', old_v_i  (ji,jj,jl)   
527                  WRITE(numout,*) ' d_v_i_trp  : ', d_v_i_trp(ji,jj,jl) , ' d_v_i_thd  : ', d_v_i_thd(ji,jj,jl) 
528                  WRITE(numout,*) ' '
529               END DO
530               inb_alp(ialert_id) = inb_alp(ialert_id) + 1
531            ENDIF
532         END DO
533      END DO
534
535      ! Alert if qns very big
536      ialert_id = 8 ! reference number of this alert
537      cl_alname(ialert_id) = ' fnsolar very big             ' ! name of the alert
538      DO jj = 1, jpj
539         DO ji = 1, jpi
540            IF(   ABS( qns(ji,jj) ) > 1500._wp  .AND.  at_i(ji,jj) > 0._wp  )  THEN
541               !
542               WRITE(numout,*) ' ALERTE 8 :   Very high non-solar heat flux'
543               WRITE(numout,*) ' ji, jj    : ', ji, jj
544               WRITE(numout,*) ' qns       : ', qns(ji,jj)
545               WRITE(numout,*) ' sst       : ', sst_m(ji,jj)
546               WRITE(numout,*) ' sss       : ', sss_m(ji,jj)
547               WRITE(numout,*) ' qcmif     : ', qcmif(ji,jj)
548               WRITE(numout,*) ' qldif     : ', qldif(ji,jj)
549               WRITE(numout,*) ' qcmif     : ', qcmif(ji,jj) / rdt_ice
550               WRITE(numout,*) ' qldif     : ', qldif(ji,jj) / rdt_ice
551               WRITE(numout,*) ' qfvbq     : ', qfvbq(ji,jj)
552               WRITE(numout,*) ' qdtcn     : ', qdtcn(ji,jj)
553               WRITE(numout,*) ' qfvbq / dt: ', qfvbq(ji,jj) / rdt_ice
554               WRITE(numout,*) ' qdtcn / dt: ', qdtcn(ji,jj) / rdt_ice
555               WRITE(numout,*) ' fdtcn     : ', fdtcn(ji,jj) 
556               WRITE(numout,*) ' fhmec     : ', fhmec(ji,jj) 
557               WRITE(numout,*) ' fheat_mec : ', fheat_mec(ji,jj) 
558               WRITE(numout,*) ' fheat_res : ', fheat_res(ji,jj) 
559               WRITE(numout,*) ' fhbri     : ', fhbri(ji,jj) 
560               !
561               CALL lim_prt_state( ji, jj, 2, '   ')
562               inb_alp(ialert_id) = inb_alp(ialert_id) + 1
563               !
564            ENDIF
565         END DO
566      END DO
567      !+++++
568 
569      ! Alert if very warm ice
570      ialert_id = 10 ! reference number of this alert
571      cl_alname(ialert_id) = ' Very warm ice                ' ! name of the alert
572      inb_alp(ialert_id) = 0
573      DO jl = 1, jpl
574         DO jk = 1, nlay_i
575            DO jj = 1, jpj
576               DO ji = 1, jpi
577                  ztmelts    =  -tmut * s_i(ji,jj,jk,jl) + rtt
578                  IF( t_i(ji,jj,jk,jl) >= ztmelts  .AND.  v_i(ji,jj,jl) > 1.e-6   &
579                     &                             .AND.  a_i(ji,jj,jl) > 0._wp   ) THEN
580                     WRITE(numout,*) ' ALERTE 10 :   Very warm ice'
581                     WRITE(numout,*) ' ji, jj, jk, jl : ', ji, jj, jk, jl
582                     WRITE(numout,*) ' t_i : ', t_i(ji,jj,jk,jl)
583                     WRITE(numout,*) ' e_i : ', e_i(ji,jj,jk,jl)
584                     WRITE(numout,*) ' s_i : ', s_i(ji,jj,jk,jl)
585                     WRITE(numout,*) ' ztmelts : ', ztmelts
586                     inb_alp(ialert_id) = inb_alp(ialert_id) + 1
587                  ENDIF
588               END DO
589            END DO
590         END DO
591      END DO
592
593      ialert_id = 1                                 ! reference number of this alert
594      cl_alname(ialert_id) = ' NO alerte 1      '   ! name of the alert
595      WRITE(numout,*)
596      WRITE(numout,*) ' All alerts at the end of ice model '
597      DO ialert_id = 1, inb_altests
598         WRITE(numout,*) ialert_id, cl_alname(ialert_id)//' : ', inb_alp(ialert_id), ' times ! '
599      END DO
600      !
601   END SUBROUTINE lim_ctl
602 
603   
604   SUBROUTINE lim_prt_state( ki, kj, kn, cd1 )
605      !!-----------------------------------------------------------------------
606      !!                   ***  ROUTINE lim_prt_state ***
607      !!                 
608      !! ** Purpose :   Writes global ice state on the (i,j) point
609      !!                in ocean.ouput
610      !!                3 possibilities exist
611      !!                n = 1/-1 -> simple ice state (plus Mechanical Check if -1)
612      !!                n = 2    -> exhaustive state
613      !!                n = 3    -> ice/ocean salt fluxes
614      !!
615      !! ** input   :   point coordinates (i,j)
616      !!                n : number of the option
617      !!-------------------------------------------------------------------
618      INTEGER         , INTENT(in) ::   ki, kj, kn    ! ocean gridpoint indices
619      CHARACTER(len=*), INTENT(in) ::   cd1           !
620      !!
621      INTEGER :: jl
622      !!-------------------------------------------------------------------
623
624      WRITE(numout,*) cd1             ! print title
625
626      !----------------
627      !  Simple state
628      !----------------
629
630      IF ( kn == 1 .OR. kn == -1 ) THEN
631         WRITE(numout,*) ' lim_prt_state - Point : ',ki,kj
632         WRITE(numout,*) ' ~~~~~~~~~~~~~~ '
633         WRITE(numout,*) ' Simple state '
634         WRITE(numout,*) ' masks s,u,v   : ', tms(ki,kj), tmu(ki,kj), tmv(ki,kj)
635         WRITE(numout,*) ' lat - long    : ', gphit(ki,kj), glamt(ki,kj)
636         WRITE(numout,*) ' Time step     : ', numit
637         WRITE(numout,*) ' - Ice drift   '
638         WRITE(numout,*) '   ~~~~~~~~~~~ '
639         WRITE(numout,*) ' u_ice(i-1,j)  : ', u_ice(ki-1,kj)
640         WRITE(numout,*) ' u_ice(i  ,j)  : ', u_ice(ki,kj)
641         WRITE(numout,*) ' v_ice(i  ,j-1): ', v_ice(ki,kj-1)
642         WRITE(numout,*) ' v_ice(i  ,j)  : ', v_ice(ki,kj)
643         WRITE(numout,*) ' strength      : ', strength(ki,kj)
644         WRITE(numout,*)
645         WRITE(numout,*) ' - Cell values '
646         WRITE(numout,*) '   ~~~~~~~~~~~ '
647         WRITE(numout,*) ' cell area     : ', area(ki,kj)
648         WRITE(numout,*) ' at_i          : ', at_i(ki,kj)       
649         WRITE(numout,*) ' vt_i          : ', vt_i(ki,kj)       
650         WRITE(numout,*) ' vt_s          : ', vt_s(ki,kj)       
651         DO jl = 1, jpl
652            WRITE(numout,*) ' - Category (', jl,')'
653            WRITE(numout,*) ' a_i           : ', a_i(ki,kj,jl)
654            WRITE(numout,*) ' ht_i          : ', ht_i(ki,kj,jl)
655            WRITE(numout,*) ' ht_s          : ', ht_s(ki,kj,jl)
656            WRITE(numout,*) ' v_i           : ', v_i(ki,kj,jl)
657            WRITE(numout,*) ' v_s           : ', v_s(ki,kj,jl)
658            WRITE(numout,*) ' e_s           : ', e_s(ki,kj,1,jl)/1.0e9
659            WRITE(numout,*) ' e_i           : ', e_i(ki,kj,1:nlay_i,jl)/1.0e9
660            WRITE(numout,*) ' t_su          : ', t_su(ki,kj,jl)
661            WRITE(numout,*) ' t_snow        : ', t_s(ki,kj,1,jl)
662            WRITE(numout,*) ' t_i           : ', t_i(ki,kj,1:nlay_i,jl)
663            WRITE(numout,*) ' sm_i          : ', sm_i(ki,kj,jl)
664            WRITE(numout,*) ' smv_i         : ', smv_i(ki,kj,jl)
665            WRITE(numout,*)
666            WRITE(numout,*) ' Pathological case : ', patho_case(ki,kj,jl)
667         END DO
668      ENDIF
669      IF( kn == -1 ) THEN
670         WRITE(numout,*) ' Mechanical Check ************** '
671         WRITE(numout,*) ' Check what means ice divergence '
672         WRITE(numout,*) ' Total ice concentration ', at_i (ki,kj)
673         WRITE(numout,*) ' Total lead fraction     ', ato_i(ki,kj)
674         WRITE(numout,*) ' Sum of both             ', ato_i(ki,kj) + at_i(ki,kj)
675         WRITE(numout,*) ' Sum of both minus 1     ', ato_i(ki,kj) + at_i(ki,kj) - 1.00
676      ENDIF
677
678
679      !--------------------
680      !  Exhaustive state
681      !--------------------
682
683      IF ( kn .EQ. 2 ) THEN
684         WRITE(numout,*) ' lim_prt_state - Point : ',ki,kj
685         WRITE(numout,*) ' ~~~~~~~~~~~~~~ '
686         WRITE(numout,*) ' Exhaustive state '
687         WRITE(numout,*) ' lat - long ', gphit(ki,kj), glamt(ki,kj)
688         WRITE(numout,*) ' Time step ', numit
689         WRITE(numout,*) 
690         WRITE(numout,*) ' - Cell values '
691         WRITE(numout,*) '   ~~~~~~~~~~~ '
692         WRITE(numout,*) ' cell area     : ', area(ki,kj)
693         WRITE(numout,*) ' at_i          : ', at_i(ki,kj)       
694         WRITE(numout,*) ' vt_i          : ', vt_i(ki,kj)       
695         WRITE(numout,*) ' vt_s          : ', vt_s(ki,kj)       
696         WRITE(numout,*) ' u_ice(i-1,j)  : ', u_ice(ki-1,kj)
697         WRITE(numout,*) ' u_ice(i  ,j)  : ', u_ice(ki,kj)
698         WRITE(numout,*) ' v_ice(i  ,j-1): ', v_ice(ki,kj-1)
699         WRITE(numout,*) ' v_ice(i  ,j)  : ', v_ice(ki,kj)
700         WRITE(numout,*) ' strength      : ', strength(ki,kj)
701         WRITE(numout,*) ' d_u_ice_dyn   : ', d_u_ice_dyn(ki,kj), ' d_v_ice_dyn   : ', d_v_ice_dyn(ki,kj)
702         WRITE(numout,*) ' old_u_ice     : ', old_u_ice(ki,kj)  , ' old_v_ice     : ', old_v_ice(ki,kj) 
703         WRITE(numout,*)
704
705         DO jl = 1, jpl
706              WRITE(numout,*) ' - Category (',jl,')'
707              WRITE(numout,*) '   ~~~~~~~~         ' 
708              WRITE(numout,*) ' ht_i       : ', ht_i(ki,kj,jl)             , ' ht_s       : ', ht_s(ki,kj,jl)
709              WRITE(numout,*) ' t_i        : ', t_i(ki,kj,1:nlay_i,jl)
710              WRITE(numout,*) ' t_su       : ', t_su(ki,kj,jl)             , ' t_s        : ', t_s(ki,kj,1,jl)
711              WRITE(numout,*) ' sm_i       : ', sm_i(ki,kj,jl)             , ' o_i        : ', o_i(ki,kj,jl)
712              WRITE(numout,*) ' a_i        : ', a_i(ki,kj,jl)              , ' old_a_i    : ', old_a_i(ki,kj,jl)   
713              WRITE(numout,*) ' d_a_i_trp  : ', d_a_i_trp(ki,kj,jl)        , ' d_a_i_thd  : ', d_a_i_thd(ki,kj,jl) 
714              WRITE(numout,*) ' v_i        : ', v_i(ki,kj,jl)              , ' old_v_i    : ', old_v_i(ki,kj,jl)   
715              WRITE(numout,*) ' d_v_i_trp  : ', d_v_i_trp(ki,kj,jl)        , ' d_v_i_thd  : ', d_v_i_thd(ki,kj,jl) 
716              WRITE(numout,*) ' v_s        : ', v_s(ki,kj,jl)              , ' old_v_s    : ', old_v_s(ki,kj,jl) 
717              WRITE(numout,*) ' d_v_s_trp  : ', d_v_s_trp(ki,kj,jl)        , ' d_v_s_thd  : ', d_v_s_thd(ki,kj,jl)
718              WRITE(numout,*) ' e_i1       : ', e_i(ki,kj,1,jl)/1.0e9      , ' old_ei1    : ', old_e_i(ki,kj,1,jl)/1.0e9 
719              WRITE(numout,*) ' de_i1_trp  : ', d_e_i_trp(ki,kj,1,jl)/1.0e9, ' de_i1_thd  : ', d_e_i_thd(ki,kj,1,jl)/1.0e9
720              WRITE(numout,*) ' e_i2       : ', e_i(ki,kj,2,jl)/1.0e9      , ' old_ei2    : ', old_e_i(ki,kj,2,jl)/1.0e9 
721              WRITE(numout,*) ' de_i2_trp  : ', d_e_i_trp(ki,kj,2,jl)/1.0e9, ' de_i2_thd  : ', d_e_i_thd(ki,kj,2,jl)/1.0e9
722              WRITE(numout,*) ' e_snow     : ', e_s(ki,kj,1,jl)            , ' old_e_snow : ', old_e_s(ki,kj,1,jl) 
723              WRITE(numout,*) ' d_e_s_trp  : ', d_e_s_trp(ki,kj,1,jl)      , ' d_e_s_thd  : ', d_e_s_thd(ki,kj,1,jl)
724              WRITE(numout,*) ' smv_i      : ', smv_i(ki,kj,jl)            , ' old_smv_i  : ', old_smv_i(ki,kj,jl)   
725              WRITE(numout,*) ' d_smv_i_trp: ', d_smv_i_trp(ki,kj,jl)      , ' d_smv_i_thd: ', d_smv_i_thd(ki,kj,jl) 
726              WRITE(numout,*) ' oa_i       : ', oa_i(ki,kj,jl)             , ' old_oa_i   : ', old_oa_i(ki,kj,jl)
727              WRITE(numout,*) ' d_oa_i_trp : ', d_oa_i_trp(ki,kj,jl)       , ' d_oa_i_thd : ', d_oa_i_thd(ki,kj,jl)
728              WRITE(numout,*) ' Path. case : ', patho_case(ki,kj,jl)
729        END DO !jl
730
731        WRITE(numout,*)
732        WRITE(numout,*) ' - Heat / FW fluxes '
733        WRITE(numout,*) '   ~~~~~~~~~~~~~~~~ '
734!       WRITE(numout,*) ' sfx_bri    : ', sfx_bri  (ki,kj)
735!       WRITE(numout,*) ' sfx        : ', sfx      (ki,kj)
736!       WRITE(numout,*) ' fsalt_res  : ', fsalt_res(ki,kj)
737        WRITE(numout,*) ' fmmec      : ', fmmec    (ki,kj)
738        WRITE(numout,*) ' fhmec      : ', fhmec    (ki,kj)
739        WRITE(numout,*) ' fhbri      : ', fhbri    (ki,kj)
740        WRITE(numout,*) ' fheat_mec  : ', fheat_mec(ki,kj)
741        WRITE(numout,*) 
742        WRITE(numout,*) ' sst        : ', sst_m(ki,kj) 
743        WRITE(numout,*) ' sss        : ', sss_m(ki,kj) 
744        WRITE(numout,*) 
745        WRITE(numout,*) ' - Stresses '
746        WRITE(numout,*) '   ~~~~~~~~ '
747        WRITE(numout,*) ' utau_ice   : ', utau_ice(ki,kj) 
748        WRITE(numout,*) ' vtau_ice   : ', vtau_ice(ki,kj)
749        WRITE(numout,*) ' utau       : ', utau    (ki,kj) 
750        WRITE(numout,*) ' vtau       : ', vtau    (ki,kj)
751        WRITE(numout,*) ' oc. vel. u : ', u_oce   (ki,kj)
752        WRITE(numout,*) ' oc. vel. v : ', v_oce   (ki,kj)
753     ENDIF
754
755     !---------------------
756     ! Salt / heat fluxes
757     !---------------------
758
759     IF ( kn .EQ. 3 ) THEN
760        WRITE(numout,*) ' lim_prt_state - Point : ',ki,kj
761        WRITE(numout,*) ' ~~~~~~~~~~~~~~ '
762        WRITE(numout,*) ' - Salt / Heat Fluxes '
763        WRITE(numout,*) '   ~~~~~~~~~~~~~~~~ '
764        WRITE(numout,*) ' lat - long ', gphit(ki,kj), glamt(ki,kj)
765        WRITE(numout,*) ' Time step ', numit
766        WRITE(numout,*)
767        WRITE(numout,*) ' - Heat fluxes at bottom interface ***'
768        WRITE(numout,*) ' qsr       : ', qsr(ki,kj)
769        WRITE(numout,*) ' qns       : ', qns(ki,kj)
770        WRITE(numout,*)
771        WRITE(numout,*) ' - Salt fluxes at bottom interface ***'
772        WRITE(numout,*) ' emp       : ', emp    (ki,kj)
773        WRITE(numout,*) ' sfx_bri   : ', sfx_bri(ki,kj)
774        WRITE(numout,*) ' sfx       : ', sfx    (ki,kj)
775        WRITE(numout,*) ' sfx_res   : ', sfx_res(ki,kj)
776        WRITE(numout,*) ' sfx_mec   : ', sfx_mec(ki,kj)
777        WRITE(numout,*) ' - Heat fluxes at bottom interface ***'
778        WRITE(numout,*) ' fheat_res : ', fheat_res(ki,kj)
779        WRITE(numout,*)
780        WRITE(numout,*) ' - Momentum fluxes '
781        WRITE(numout,*) ' utau      : ', utau(ki,kj) 
782        WRITE(numout,*) ' vtau      : ', vtau(ki,kj)
783      ENDIF
784      WRITE(numout,*) ' '
785      !
786   END SUBROUTINE lim_prt_state
787
788#else
789   !!----------------------------------------------------------------------
790   !!   Default option           Dummy module      NO LIM 3.0 sea-ice model
791   !!----------------------------------------------------------------------
792CONTAINS
793   SUBROUTINE sbc_ice_lim ( kt, kblk )     ! Dummy routine
794      WRITE(*,*) 'sbc_ice_lim: You should not have seen this print! error?', kt, kblk
795   END SUBROUTINE sbc_ice_lim
796#endif
797
798   !!======================================================================
799END MODULE sbcice_lim
Note: See TracBrowser for help on using the repository browser.