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

source: branches/dev_r4028_CNRS_LIM3_MV2014/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90 @ 4358

Last change on this file since 4358 was 4358, checked in by vancop, 10 years ago

[test enhanced CORE over sea ice, see ticket #1219]

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