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_r4650_UKMO13_CICE_changes_take2/NEMOGCM/NEMO/OPA_SRC/SBC – NEMO

source: branches/2014/dev_r4650_UKMO13_CICE_changes_take2/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90 @ 4921

Last change on this file since 4921 was 4921, checked in by timgraham, 9 years ago

merged with revision 4879 of trunk

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