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.
Changeset 6092 for branches/2015/dev_agrif_v3_6_STABLE/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90 – NEMO

Ignore:
Timestamp:
2015-12-17T15:19:01+01:00 (9 years ago)
Author:
timgraham
Message:

Merged in trunk at r5518 (branch point of 3.6 stable)

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2015/dev_agrif_v3_6_STABLE/NEMOGCM/NEMO/OPA_SRC/SBC/sbcice_lim.F90

    r4765 r6092  
    1212   !!            3.4  ! 2011-01  (A Porter)  dynamical allocation 
    1313   !!             -   ! 2012-10  (C. Rousset)  add lim_diahsb 
     14   !!            3.6  ! 2014-07  (M. Vancoppenolle, G. Madec, O. Marti) revise coupled interface 
    1415   !!---------------------------------------------------------------------- 
    1516#if defined key_lim3 
     
    1819   !!---------------------------------------------------------------------- 
    1920   !!   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 
    2221   !!---------------------------------------------------------------------- 
    2322   USE oce             ! ocean dynamics and tracers 
    2423   USE dom_oce         ! ocean space and time domain 
    25    USE par_ice         ! sea-ice parameters 
    2624   USE ice             ! LIM-3: ice variables 
    27    USE iceini          ! LIM-3: ice initialisation 
     25   USE thd_ice         ! LIM-3: thermodynamical variables 
    2826   USE dom_ice         ! LIM-3: ice domain 
    2927 
     
    3937   USE limdyn          ! Ice dynamics 
    4038   USE limtrp          ! Ice transport 
     39   USE limhdf          ! Ice horizontal diffusion 
    4140   USE limthd          ! Ice thermodynamics 
    42    USE limitd_th       ! Thermodynamics on ice thickness distribution  
    4341   USE limitd_me       ! Mechanics on ice thickness distribution 
    4442   USE limsbc          ! sea surface boundary condition 
     
    4644   USE limwri          ! Ice outputs 
    4745   USE limrst          ! Ice restarts 
    48    USE limupdate1       ! update of global variables 
    49    USE limupdate2       ! update of global variables 
     46   USE limupdate1      ! update of global variables 
     47   USE limupdate2      ! update of global variables 
    5048   USE limvar          ! Ice variables switch 
     49 
     50   USE limmsh          ! LIM mesh 
     51   USE limistate       ! LIM initial state 
     52   USE limthd_sal      ! LIM ice thermodynamics: salinity 
    5153 
    5254   USE c1d             ! 1D vertical configuration 
     
    5961   USE prtctl          ! Print control 
    6062   USE lib_fortran     !  
    61    USE cpl_oasis3, ONLY : lk_cpl 
     63   USE limctl 
    6264 
    6365#if defined key_bdy  
     
    6971 
    7072   PUBLIC sbc_ice_lim  ! routine called by sbcmod.F90 
    71    PUBLIC lim_prt_state 
     73   PUBLIC sbc_lim_init ! routine called by sbcmod.F90 
    7274    
    7375   !! * Substitutions 
     
    8183CONTAINS 
    8284 
    83    FUNCTION fice_cell_ave ( ptab) 
     85   !!====================================================================== 
     86 
     87   SUBROUTINE sbc_ice_lim( kt, kblk ) 
     88      !!--------------------------------------------------------------------- 
     89      !!                  ***  ROUTINE sbc_ice_lim  *** 
     90      !!                    
     91      !! ** Purpose :   update the ocean surface boundary condition via the  
     92      !!                Louvain la Neuve Sea Ice Model time stepping  
     93      !! 
     94      !! ** Method  :   ice model time stepping 
     95      !!              - call the ice dynamics routine  
     96      !!              - call the ice advection/diffusion routine  
     97      !!              - call the ice thermodynamics routine  
     98      !!              - call the routine that computes mass and  
     99      !!                heat fluxes at the ice/ocean interface 
     100      !!              - save the outputs  
     101      !!              - save the outputs for restart when necessary 
     102      !! 
     103      !! ** Action  : - time evolution of the LIM sea-ice model 
     104      !!              - update all sbc variables below sea-ice: 
     105      !!                utau, vtau, taum, wndm, qns , qsr, emp , sfx  
     106      !!--------------------------------------------------------------------- 
     107      INTEGER, INTENT(in) ::   kt      ! ocean time step 
     108      INTEGER, INTENT(in) ::   kblk    ! type of bulk (=3 CLIO, =4 CORE, =5 COUPLED) 
     109      !! 
     110      INTEGER  ::   jl                 ! dummy loop index 
     111      REAL(wp), POINTER, DIMENSION(:,:,:)   ::   zalb_os, zalb_cs  ! ice albedo under overcast/clear sky 
     112      REAL(wp), POINTER, DIMENSION(:,:,:)   ::   zalb_ice          ! mean ice albedo (for coupled) 
     113      REAL(wp), POINTER, DIMENSION(:,:  )   ::   zutau_ice, zvtau_ice  
     114      !!---------------------------------------------------------------------- 
     115 
     116      IF( nn_timing == 1 )  CALL timing_start('sbc_ice_lim') 
     117 
     118      !-----------------------! 
     119      ! --- Ice time step --- ! 
     120      !-----------------------! 
     121      IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN 
     122 
     123         ! mean surface ocean current at ice velocity point (C-grid dynamics :  U- & V-points as the ocean) 
     124         u_oce(:,:) = ssu_m(:,:) * umask(:,:,1) 
     125         v_oce(:,:) = ssv_m(:,:) * vmask(:,:,1) 
     126          
     127         ! masked sea surface freezing temperature [Kelvin] (set to rt0 over land) 
     128         t_bo(:,:) = ( eos_fzp( sss_m ) + rt0 ) * tmask(:,:,1) + rt0 * ( 1._wp - tmask(:,:,1) )   
     129          
     130         ! Mask sea ice surface temperature (set to rt0 over land) 
     131         DO jl = 1, jpl 
     132            t_su(:,:,jl) = t_su(:,:,jl) * tmask(:,:,1) + rt0 * ( 1._wp - tmask(:,:,1) ) 
     133         END DO      
     134         ! 
     135         !------------------------------------------------!                                            
     136         ! --- Dynamical coupling with the atmosphere --- !                                            
     137         !------------------------------------------------! 
     138         ! It provides the following fields: 
     139         ! utau_ice, vtau_ice : surface ice stress (U- & V-points)   [N/m2] 
     140         !----------------------------------------------------------------- 
     141         SELECT CASE( kblk ) 
     142         CASE( jp_clio    )   ;   CALL blk_ice_clio_tau                         ! CLIO bulk formulation             
     143         CASE( jp_core    )   ;   CALL blk_ice_core_tau                         ! CORE bulk formulation 
     144         CASE( jp_purecpl )   ;   CALL sbc_cpl_ice_tau( utau_ice , vtau_ice )   ! Coupled   formulation 
     145         END SELECT 
     146          
     147         IF( ln_mixcpl) THEN   ! Case of a mixed Bulk/Coupled formulation 
     148            CALL wrk_alloc( jpi,jpj    , zutau_ice, zvtau_ice) 
     149            CALL sbc_cpl_ice_tau( zutau_ice , zvtau_ice ) 
     150            utau_ice(:,:) = utau_ice(:,:) * xcplmask(:,:,0) + zutau_ice(:,:) * ( 1. - xcplmask(:,:,0) ) 
     151            vtau_ice(:,:) = vtau_ice(:,:) * xcplmask(:,:,0) + zvtau_ice(:,:) * ( 1. - xcplmask(:,:,0) ) 
     152            CALL wrk_dealloc( jpi,jpj  , zutau_ice, zvtau_ice) 
     153         ENDIF 
     154 
     155         !-------------------------------------------------------! 
     156         ! --- ice dynamics and transport (except in 1D case) ---! 
     157         !-------------------------------------------------------! 
     158         numit = numit + nn_fsbc                  ! Ice model time step 
     159         !                                                    
     160         CALL sbc_lim_bef                         ! Store previous ice values 
     161         CALL sbc_lim_diag0                       ! set diag of mass, heat and salt fluxes to 0 
     162         CALL lim_rst_opn( kt )                   ! Open Ice restart file 
     163         ! 
     164         IF( .NOT. lk_c1d ) THEN 
     165            ! 
     166            CALL lim_dyn( kt )                    ! Ice dynamics    ( rheology/dynamics )    
     167            ! 
     168            CALL lim_trp( kt )                    ! Ice transport   ( Advection/diffusion ) 
     169            ! 
     170            IF( nn_monocat /= 2 ) CALL lim_itd_me ! Mechanical redistribution ! (ridging/rafting) 
     171            ! 
     172#if defined key_bdy 
     173            CALL bdy_ice_lim( kt )                ! bdy ice thermo  
     174            IF( ln_icectl )       CALL lim_prt( kt, iiceprt, jiceprt, 1, ' - ice thermo bdy - ' ) 
     175#endif 
     176            ! 
     177            CALL lim_update1( kt )                ! Corrections 
     178            ! 
     179         ENDIF 
     180          
     181         ! previous lead fraction and ice volume for flux calculations 
     182         CALL sbc_lim_bef                         
     183         CALL lim_var_glo2eqv                     ! ht_i and ht_s for ice albedo calculation 
     184         CALL lim_var_agg(1)                      ! at_i for coupling (via pfrld)  
     185         pfrld(:,:)   = 1._wp - at_i(:,:) 
     186         phicif(:,:)  = vt_i(:,:) 
     187          
     188         !------------------------------------------------------!                                            
     189         ! --- Thermodynamical coupling with the atmosphere --- !                                            
     190         !------------------------------------------------------! 
     191         ! It provides the following fields: 
     192         ! qsr_ice , qns_ice  : solar & non solar heat flux over ice   (T-point)         [W/m2] 
     193         ! qla_ice            : latent heat flux over ice              (T-point)         [W/m2] 
     194         ! dqns_ice, dqla_ice : non solar & latent heat sensistivity   (T-point)         [W/m2] 
     195         ! tprecip , sprecip  : total & solid precipitation            (T-point)         [Kg/m2/s] 
     196         ! fr1_i0  , fr2_i0   : 1sr & 2nd fraction of qsr penetration in ice             [%] 
     197         !---------------------------------------------------------------------------------------- 
     198         CALL wrk_alloc( jpi,jpj,jpl, zalb_os, zalb_cs, zalb_ice ) 
     199         CALL albedo_ice( t_su, ht_i, ht_s, zalb_cs, zalb_os ) ! cloud-sky and overcast-sky ice albedos 
     200 
     201         SELECT CASE( kblk ) 
     202         CASE( jp_clio )                                       ! CLIO bulk formulation 
     203            ! In CLIO the cloud fraction is read in the climatology and the all-sky albedo  
     204            ! (zalb_ice) is computed within the bulk routine 
     205            CALL blk_ice_clio_flx( t_su, zalb_cs, zalb_os, zalb_ice ) 
     206            IF( ln_mixcpl      ) CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=zalb_ice, psst=sst_m, pist=t_su ) 
     207            IF( nn_limflx /= 2 ) CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx ) 
     208         CASE( jp_core )                                       ! CORE bulk formulation 
     209            ! albedo depends on cloud fraction because of non-linear spectral effects 
     210            zalb_ice(:,:,:) = ( 1. - cldf_ice ) * zalb_cs(:,:,:) + cldf_ice * zalb_os(:,:,:) 
     211            CALL blk_ice_core_flx( t_su, zalb_ice ) 
     212            IF( ln_mixcpl      ) CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=zalb_ice, psst=sst_m, pist=t_su ) 
     213            IF( nn_limflx /= 2 ) CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx ) 
     214         CASE ( jp_purecpl ) 
     215            ! albedo depends on cloud fraction because of non-linear spectral effects 
     216            zalb_ice(:,:,:) = ( 1. - cldf_ice ) * zalb_cs(:,:,:) + cldf_ice * zalb_os(:,:,:) 
     217                                 CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=zalb_ice, psst=sst_m, pist=t_su ) 
     218            ! clem: evap_ice is forced to 0 in coupled mode for now  
     219            !       but it needs to be changed (along with modif in limthd_dh) once heat flux from evap will be avail. from atm. models 
     220            evap_ice  (:,:,:) = 0._wp   ;   devap_ice (:,:,:) = 0._wp 
     221            IF( nn_limflx == 2 ) CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice, dqns_ice, evap_ice, devap_ice, nn_limflx ) 
     222         END SELECT 
     223         CALL wrk_dealloc( jpi,jpj,jpl, zalb_os, zalb_cs, zalb_ice ) 
     224 
     225         !----------------------------! 
     226         ! --- ice thermodynamics --- ! 
     227         !----------------------------! 
     228         CALL lim_thd( kt )                         ! Ice thermodynamics       
     229         ! 
     230         CALL lim_update2( kt )                     ! Corrections 
     231         ! 
     232         CALL lim_sbc_flx( kt )                     ! Update surface ocean mass, heat and salt fluxes 
     233         ! 
     234         IF(ln_limdiaout) CALL lim_diahsb           ! Diagnostics and outputs  
     235         ! 
     236         CALL lim_wri( 1 )                          ! Ice outputs  
     237         ! 
     238         IF( kt == nit000 .AND. ln_rstart )   & 
     239            &             CALL iom_close( numrir )  ! close input ice restart file 
     240         ! 
     241         IF( lrst_ice )   CALL lim_rst_write( kt )  ! Ice restart file  
     242         ! 
     243         IF( ln_icectl )  CALL lim_ctl( kt )        ! alerts in case of model crash 
     244         ! 
     245      ENDIF   ! End sea-ice time step only 
     246 
     247      !-------------------------! 
     248      ! --- Ocean time step --- ! 
     249      !-------------------------! 
     250      ! Update surface ocean stresses (only in ice-dynamic case) otherwise the atm.-ocean stresses are used everywhere 
     251      IF( ln_limdyn )     CALL lim_sbc_tau( kt, ub(:,:,1), vb(:,:,1) )  ! using before instantaneous surf. currents 
     252!!gm   remark, the ocean-ice stress is not saved in ice diag call above .....  find a solution!!! 
     253      ! 
     254      IF( nn_timing == 1 ) CALL timing_stop('sbc_ice_lim') 
     255      ! 
     256   END SUBROUTINE sbc_ice_lim 
     257    
     258 
     259   SUBROUTINE sbc_lim_init 
     260      !!---------------------------------------------------------------------- 
     261      !!                  ***  ROUTINE sbc_lim_init  *** 
     262      !! 
     263      !! ** purpose :   Allocate all the dynamic arrays of the LIM-3 modules 
     264      !!---------------------------------------------------------------------- 
     265      INTEGER :: ierr 
     266      !!---------------------------------------------------------------------- 
     267      IF(lwp) WRITE(numout,*) 
     268      IF(lwp) WRITE(numout,*) 'sbc_ice_lim : update ocean surface boudary condition'  
     269      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~   via Louvain la Neuve Ice Model (LIM-3) time stepping' 
     270      ! 
     271                                       ! Open the reference and configuration namelist files and namelist output file  
     272      CALL ctl_opn( numnam_ice_ref, 'namelist_ice_ref',    'OLD',     'FORMATTED', 'SEQUENTIAL', -1, numout, lwp )  
     273      CALL ctl_opn( numnam_ice_cfg, 'namelist_ice_cfg',    'OLD',     'FORMATTED', 'SEQUENTIAL', -1, numout, lwp ) 
     274      IF(lwm) CALL ctl_opn( numoni, 'output.namelist.ice', 'UNKNOWN', 'FORMATTED', 'SEQUENTIAL', -1, numout, lwp, 1 ) 
     275 
     276      CALL ice_run                     ! set some ice run parameters 
     277      ! 
     278      !                                ! Allocate the ice arrays 
     279      ierr =        ice_alloc        ()      ! ice variables 
     280      ierr = ierr + dom_ice_alloc    ()      ! domain 
     281      ierr = ierr + sbc_ice_alloc    ()      ! surface forcing 
     282      ierr = ierr + thd_ice_alloc    ()      ! thermodynamics 
     283      ierr = ierr + lim_itd_me_alloc ()      ! ice thickness distribution - mechanics 
     284      ! 
     285      IF( lk_mpp    )   CALL mpp_sum( ierr ) 
     286      IF( ierr /= 0 )   CALL ctl_stop('STOP', 'sbc_lim_init : unable to allocate ice arrays') 
     287      ! 
     288      !                                ! adequation jpk versus ice/snow layers/categories 
     289      IF( jpl > jpk .OR. (nlay_i+1) > jpk .OR. nlay_s > jpk )   & 
     290         &      CALL ctl_stop( 'STOP',                          & 
     291         &     'sbc_lim_init: the 3rd dimension of workspace arrays is too small.',   & 
     292         &     'use more ocean levels or less ice/snow layers/categories.' ) 
     293      ! 
     294      CALL lim_itd_init                ! ice thickness distribution initialization 
     295      ! 
     296      CALL lim_hdf_init                ! set ice horizontal diffusion computation parameters 
     297      ! 
     298      CALL lim_thd_init                ! set ice thermodynics parameters 
     299      ! 
     300      CALL lim_thd_sal_init            ! set ice salinity parameters 
     301      ! 
     302      CALL lim_msh                     ! ice mesh initialization 
     303      ! 
     304      CALL lim_itd_me_init             ! ice thickness distribution initialization for mecanical deformation 
     305      !                                ! Initial sea-ice state 
     306      IF( .NOT. ln_rstart ) THEN              ! start from rest: sea-ice deduced from sst 
     307         numit = 0 
     308         numit = nit000 - 1 
     309         CALL lim_istate 
     310      ELSE                                    ! start from a restart file 
     311         CALL lim_rst_read 
     312         numit = nit000 - 1 
     313      ENDIF 
     314      CALL lim_var_agg(1) 
     315      CALL lim_var_glo2eqv 
     316      ! 
     317      CALL lim_sbc_init                 ! ice surface boundary condition    
     318      ! 
     319      fr_i(:,:)     = at_i(:,:)         ! initialisation of sea-ice fraction 
     320      tn_ice(:,:,:) = t_su(:,:,:)       ! initialisation of surface temp for coupled simu 
     321      ! 
     322      nstart = numit  + nn_fsbc       
     323      nitrun = nitend - nit000 + 1  
     324      nlast  = numit  + nitrun  
     325      ! 
     326      IF( nstock == 0 )   nstock = nlast + 1 
     327      ! 
     328   END SUBROUTINE sbc_lim_init 
     329 
     330 
     331   SUBROUTINE ice_run 
     332      !!------------------------------------------------------------------- 
     333      !!                  ***  ROUTINE ice_run *** 
     334      !!                  
     335      !! ** Purpose :   Definition some run parameter for ice model 
     336      !! 
     337      !! ** Method  :   Read the namicerun namelist and check the parameter  
     338      !!              values called at the first timestep (nit000) 
     339      !! 
     340      !! ** input   :   Namelist namicerun 
     341      !!------------------------------------------------------------------- 
     342      INTEGER  ::   ios                 ! Local integer output status for namelist read 
     343      NAMELIST/namicerun/ jpl, nlay_i, nlay_s, cn_icerst_in, cn_icerst_indir, cn_icerst_out, cn_icerst_outdir,  & 
     344         &                ln_limdyn, rn_amax, ln_limdiahsb, ln_limdiaout, ln_icectl, iiceprt, jiceprt   
     345      !!------------------------------------------------------------------- 
     346      !                     
     347      REWIND( numnam_ice_ref )              ! Namelist namicerun in reference namelist : Parameters for ice 
     348      READ  ( numnam_ice_ref, namicerun, IOSTAT = ios, ERR = 901) 
     349901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicerun in reference namelist', lwp ) 
     350 
     351      REWIND( numnam_ice_cfg )              ! Namelist namicerun in configuration namelist : Parameters for ice 
     352      READ  ( numnam_ice_cfg, namicerun, IOSTAT = ios, ERR = 902 ) 
     353902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namicerun in configuration namelist', lwp ) 
     354      IF(lwm) WRITE ( numoni, namicerun ) 
     355      ! 
     356      ! 
     357      IF(lwp) THEN                        ! control print 
     358         WRITE(numout,*) 
     359         WRITE(numout,*) 'ice_run : ice share parameters for dynamics/advection/thermo of sea-ice' 
     360         WRITE(numout,*) ' ~~~~~~' 
     361         WRITE(numout,*) '   number of ice  categories                               = ', jpl 
     362         WRITE(numout,*) '   number of ice  layers                                   = ', nlay_i 
     363         WRITE(numout,*) '   number of snow layers                                   = ', nlay_s 
     364         WRITE(numout,*) '   switch for ice dynamics (1) or not (0)      ln_limdyn   = ', ln_limdyn 
     365         WRITE(numout,*) '   maximum ice concentration                               = ', rn_amax  
     366         WRITE(numout,*) '   Diagnose heat/salt budget or not          ln_limdiahsb  = ', ln_limdiahsb 
     367         WRITE(numout,*) '   Output   heat/salt budget or not          ln_limdiaout  = ', ln_limdiaout 
     368         WRITE(numout,*) '   control prints in ocean.out for (i,j)=(iiceprt,jiceprt) = ', ln_icectl 
     369         WRITE(numout,*) '   i-index for control prints (ln_icectl=true)             = ', iiceprt 
     370         WRITE(numout,*) '   j-index for control prints (ln_icectl=true)             = ', jiceprt 
     371      ENDIF 
     372      ! 
     373      ! sea-ice timestep and inverse 
     374      rdt_ice   = nn_fsbc * rdttra(1)   
     375      r1_rdtice = 1._wp / rdt_ice  
     376 
     377      ! inverse of nlay_i and nlay_s 
     378      r1_nlay_i = 1._wp / REAL( nlay_i, wp ) 
     379      r1_nlay_s = 1._wp / REAL( nlay_s, wp ) 
     380      ! 
     381#if defined key_bdy 
     382      IF( lwp .AND. ln_limdiahsb )  CALL ctl_warn('online conservation check activated but it does not work with BDY') 
     383#endif 
     384      ! 
     385   END SUBROUTINE ice_run 
     386 
     387 
     388   SUBROUTINE lim_itd_init 
     389      !!------------------------------------------------------------------ 
     390      !!                ***  ROUTINE lim_itd_init *** 
     391      !! 
     392      !! ** Purpose :   Initializes the ice thickness distribution 
     393      !! ** Method  :   ... 
     394      !! ** input   :   Namelist namiceitd 
     395      !!------------------------------------------------------------------- 
     396      INTEGER  ::   ios                 ! Local integer output status for namelist read 
     397      NAMELIST/namiceitd/ nn_catbnd, rn_himean 
     398      ! 
     399      INTEGER  ::   jl                   ! dummy loop index 
     400      REAL(wp) ::   zc1, zc2, zc3, zx1   ! local scalars 
     401      REAL(wp) ::   zhmax, znum, zden, zalpha ! 
     402      !!------------------------------------------------------------------ 
     403      ! 
     404      REWIND( numnam_ice_ref )              ! Namelist namiceitd in reference namelist : Parameters for ice 
     405      READ  ( numnam_ice_ref, namiceitd, IOSTAT = ios, ERR = 903) 
     406903   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namiceitd in reference namelist', lwp ) 
     407 
     408      REWIND( numnam_ice_cfg )              ! Namelist namiceitd in configuration namelist : Parameters for ice 
     409      READ  ( numnam_ice_cfg, namiceitd, IOSTAT = ios, ERR = 904 ) 
     410904   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namiceitd in configuration namelist', lwp ) 
     411      IF(lwm) WRITE ( numoni, namiceitd ) 
     412      ! 
     413      ! 
     414      IF(lwp) THEN                        ! control print 
     415         WRITE(numout,*) 
     416         WRITE(numout,*) 'ice_itd : ice cat distribution' 
     417         WRITE(numout,*) ' ~~~~~~' 
     418         WRITE(numout,*) '   shape of ice categories distribution                          nn_catbnd = ', nn_catbnd 
     419         WRITE(numout,*) '   mean ice thickness in the domain (only active if nn_catbnd=2) rn_himean = ', rn_himean 
     420      ENDIF 
     421 
     422      !---------------------------------- 
     423      !- Thickness categories boundaries  
     424      !---------------------------------- 
     425      IF(lwp) WRITE(numout,*) 
     426      IF(lwp) WRITE(numout,*) 'lim_itd_init : Initialization of ice cat distribution ' 
     427      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~' 
     428 
     429      hi_max(:) = 0._wp 
     430 
     431      SELECT CASE ( nn_catbnd  )        
     432                                   !---------------------- 
     433         CASE (1)                  ! tanh function (CICE) 
     434                                   !---------------------- 
     435         zc1 =  3._wp / REAL( jpl, wp ) 
     436         zc2 = 10._wp * zc1 
     437         zc3 =  3._wp 
     438 
     439         DO jl = 1, jpl 
     440            zx1 = REAL( jl-1, wp ) / REAL( jpl, wp ) 
     441            hi_max(jl) = hi_max(jl-1) + zc1 + zc2 * (1._wp + TANH( zc3 * (zx1 - 1._wp ) ) ) 
     442         END DO 
     443 
     444                                   !---------------------- 
     445         CASE (2)                  ! h^(-alpha) function 
     446                                   !---------------------- 
     447         zalpha = 0.05             ! exponent of the transform function 
     448 
     449         zhmax  = 3.*rn_himean 
     450 
     451         DO jl = 1, jpl  
     452            znum = jpl * ( zhmax+1 )**zalpha 
     453            zden = ( jpl - jl ) * ( zhmax+1 )**zalpha + jl 
     454            hi_max(jl) = ( znum / zden )**(1./zalpha) - 1 
     455         END DO 
     456 
     457      END SELECT 
     458 
     459      DO jl = 1, jpl 
     460         hi_mean(jl) = ( hi_max(jl) + hi_max(jl-1) ) * 0.5_wp 
     461      END DO 
     462 
     463      ! Set hi_max(jpl) to a big value to ensure that all ice is thinner than hi_max(jpl) 
     464      hi_max(jpl) = 99._wp 
     465 
     466      IF(lwp) WRITE(numout,*) ' Thickness category boundaries ' 
     467      IF(lwp) WRITE(numout,*) ' hi_max ', hi_max(0:jpl) 
     468      ! 
     469   END SUBROUTINE lim_itd_init 
     470 
     471    
     472   SUBROUTINE ice_lim_flx( ptn_ice, palb_ice, pqns_ice, pqsr_ice, pdqn_ice, pevap_ice, pdevap_ice, k_limflx ) 
     473      !!--------------------------------------------------------------------- 
     474      !!                  ***  ROUTINE ice_lim_flx  *** 
     475      !!                    
     476      !! ** Purpose :   update the ice surface boundary condition by averaging and / or 
     477      !!                redistributing fluxes on ice categories                    
     478      !! 
     479      !! ** Method  :   average then redistribute  
     480      !! 
     481      !! ** Action  :    
     482      !!--------------------------------------------------------------------- 
     483      INTEGER                   , INTENT(in   ) ::   k_limflx   ! =-1 do nothing; =0 average ;  
     484                                                                ! =1 average and redistribute ; =2 redistribute 
     485      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   ptn_ice    ! ice surface temperature  
     486      REAL(wp), DIMENSION(:,:,:), INTENT(in   ) ::   palb_ice   ! ice albedo 
     487      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pqns_ice   ! non solar flux 
     488      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pqsr_ice   ! net solar flux 
     489      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pdqn_ice   ! non solar flux sensitivity 
     490      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pevap_ice  ! sublimation 
     491      REAL(wp), DIMENSION(:,:,:), INTENT(inout) ::   pdevap_ice ! sublimation sensitivity 
     492      ! 
     493      INTEGER  ::   jl      ! dummy loop index 
     494      ! 
     495      REAL(wp), POINTER, DIMENSION(:,:) :: zalb_m    ! Mean albedo over all categories 
     496      REAL(wp), POINTER, DIMENSION(:,:) :: ztem_m    ! Mean temperature over all categories 
     497      ! 
     498      REAL(wp), POINTER, DIMENSION(:,:) :: z_qsr_m   ! Mean solar heat flux over all categories 
     499      REAL(wp), POINTER, DIMENSION(:,:) :: z_qns_m   ! Mean non solar heat flux over all categories 
     500      REAL(wp), POINTER, DIMENSION(:,:) :: z_evap_m  ! Mean sublimation over all categories 
     501      REAL(wp), POINTER, DIMENSION(:,:) :: z_dqn_m   ! Mean d(qns)/dT over all categories 
     502      REAL(wp), POINTER, DIMENSION(:,:) :: z_devap_m ! Mean d(evap)/dT over all categories 
     503      !!---------------------------------------------------------------------- 
     504 
     505      IF( nn_timing == 1 )  CALL timing_start('ice_lim_flx') 
     506      ! 
     507      ! 
     508      SELECT CASE( k_limflx )                              !==  averaged on all ice categories  ==! 
     509      CASE( 0 , 1 ) 
     510         CALL wrk_alloc( jpi,jpj, z_qsr_m, z_qns_m, z_evap_m, z_dqn_m, z_devap_m) 
     511         ! 
     512         z_qns_m  (:,:) = fice_ice_ave ( pqns_ice (:,:,:) ) 
     513         z_qsr_m  (:,:) = fice_ice_ave ( pqsr_ice (:,:,:) ) 
     514         z_dqn_m  (:,:) = fice_ice_ave ( pdqn_ice (:,:,:) ) 
     515         z_evap_m (:,:) = fice_ice_ave ( pevap_ice (:,:,:) ) 
     516         z_devap_m(:,:) = fice_ice_ave ( pdevap_ice (:,:,:) ) 
     517         DO jl = 1, jpl 
     518            pdqn_ice  (:,:,jl) = z_dqn_m(:,:) 
     519            pdevap_ice(:,:,jl) = z_devap_m(:,:) 
     520         END DO 
     521         ! 
     522         DO jl = 1, jpl 
     523            pqns_ice (:,:,jl) = z_qns_m(:,:) 
     524            pqsr_ice (:,:,jl) = z_qsr_m(:,:) 
     525            pevap_ice(:,:,jl) = z_evap_m(:,:) 
     526         END DO 
     527         ! 
     528         CALL wrk_dealloc( jpi,jpj, z_qsr_m, z_qns_m, z_evap_m, z_dqn_m, z_devap_m) 
     529      END SELECT 
     530 
     531      SELECT CASE( k_limflx )                              !==  redistribution on all ice categories  ==! 
     532      CASE( 1 , 2 ) 
     533         CALL wrk_alloc( jpi,jpj, zalb_m, ztem_m ) 
     534         ! 
     535         zalb_m(:,:) = fice_ice_ave ( palb_ice (:,:,:) )  
     536         ztem_m(:,:) = fice_ice_ave ( ptn_ice  (:,:,:) )  
     537         DO jl = 1, jpl 
     538            pqns_ice (:,:,jl) = pqns_ice (:,:,jl) + pdqn_ice  (:,:,jl) * ( ptn_ice(:,:,jl) - ztem_m(:,:) ) 
     539            pevap_ice(:,:,jl) = pevap_ice(:,:,jl) + pdevap_ice(:,:,jl) * ( ptn_ice(:,:,jl) - ztem_m(:,:) ) 
     540            pqsr_ice (:,:,jl) = pqsr_ice (:,:,jl) * ( 1._wp - palb_ice(:,:,jl) ) / ( 1._wp - zalb_m(:,:) )  
     541         END DO 
     542         ! 
     543         CALL wrk_dealloc( jpi,jpj, zalb_m, ztem_m ) 
     544      END SELECT 
     545      ! 
     546      IF( nn_timing == 1 )  CALL timing_stop('ice_lim_flx') 
     547      ! 
     548   END SUBROUTINE ice_lim_flx 
     549 
     550   SUBROUTINE sbc_lim_bef 
     551      !!---------------------------------------------------------------------- 
     552      !!                  ***  ROUTINE sbc_lim_bef  *** 
     553      !! 
     554      !! ** purpose :  store ice variables at "before" time step  
     555      !!---------------------------------------------------------------------- 
     556      a_i_b  (:,:,:)   = a_i  (:,:,:)     ! ice area 
     557      e_i_b  (:,:,:,:) = e_i  (:,:,:,:)   ! ice thermal energy 
     558      v_i_b  (:,:,:)   = v_i  (:,:,:)     ! ice volume 
     559      v_s_b  (:,:,:)   = v_s  (:,:,:)     ! snow volume  
     560      e_s_b  (:,:,:,:) = e_s  (:,:,:,:)   ! snow thermal energy 
     561      smv_i_b(:,:,:)   = smv_i(:,:,:)     ! salt content 
     562      oa_i_b (:,:,:)   = oa_i (:,:,:)     ! areal age content 
     563      u_ice_b(:,:)     = u_ice(:,:) 
     564      v_ice_b(:,:)     = v_ice(:,:) 
     565       
     566   END SUBROUTINE sbc_lim_bef 
     567 
     568   SUBROUTINE sbc_lim_diag0 
     569      !!---------------------------------------------------------------------- 
     570      !!                  ***  ROUTINE sbc_lim_diag0  *** 
     571      !! 
     572      !! ** purpose :  set ice-ocean and ice-atm. fluxes to zeros at the beggining 
     573      !!               of the time step 
     574      !!---------------------------------------------------------------------- 
     575      sfx    (:,:) = 0._wp   ; 
     576      sfx_bri(:,:) = 0._wp   ;  
     577      sfx_sni(:,:) = 0._wp   ;   sfx_opw(:,:) = 0._wp 
     578      sfx_bog(:,:) = 0._wp   ;   sfx_dyn(:,:) = 0._wp 
     579      sfx_bom(:,:) = 0._wp   ;   sfx_sum(:,:) = 0._wp 
     580      sfx_res(:,:) = 0._wp 
     581       
     582      wfx_snw(:,:) = 0._wp   ;   wfx_ice(:,:) = 0._wp 
     583      wfx_sni(:,:) = 0._wp   ;   wfx_opw(:,:) = 0._wp 
     584      wfx_bog(:,:) = 0._wp   ;   wfx_dyn(:,:) = 0._wp 
     585      wfx_bom(:,:) = 0._wp   ;   wfx_sum(:,:) = 0._wp 
     586      wfx_res(:,:) = 0._wp   ;   wfx_sub(:,:) = 0._wp 
     587      wfx_spr(:,:) = 0._wp   ;    
     588       
     589      hfx_thd(:,:) = 0._wp   ;    
     590      hfx_snw(:,:) = 0._wp   ;   hfx_opw(:,:) = 0._wp 
     591      hfx_bog(:,:) = 0._wp   ;   hfx_dyn(:,:) = 0._wp 
     592      hfx_bom(:,:) = 0._wp   ;   hfx_sum(:,:) = 0._wp 
     593      hfx_res(:,:) = 0._wp   ;   hfx_sub(:,:) = 0._wp 
     594      hfx_spr(:,:) = 0._wp   ;   hfx_dif(:,:) = 0._wp  
     595      hfx_err(:,:) = 0._wp   ;   hfx_err_rem(:,:) = 0._wp 
     596      hfx_err_dif(:,:) = 0._wp   ; 
     597 
     598      afx_tot(:,:) = 0._wp   ; 
     599      afx_dyn(:,:) = 0._wp   ;   afx_thd(:,:) = 0._wp 
     600 
     601      diag_heat(:,:) = 0._wp ;   diag_smvi(:,:) = 0._wp ; 
     602      diag_vice(:,:) = 0._wp ;   diag_vsnw(:,:) = 0._wp ; 
     603       
     604   END SUBROUTINE sbc_lim_diag0 
     605 
     606      
     607   FUNCTION fice_cell_ave ( ptab ) 
    84608      !!-------------------------------------------------------------------------- 
    85609      !! * Compute average over categories, for grid cell (ice covered and free ocean) 
     
    90614       
    91615      fice_cell_ave (:,:) = 0.0_wp 
    92        
    93616      DO jl = 1, jpl 
    94          fice_cell_ave (:,:) = fice_cell_ave (:,:) & 
    95             &                  + a_i (:,:,jl) * ptab (:,:,jl) 
     617         fice_cell_ave (:,:) = fice_cell_ave (:,:) + a_i (:,:,jl) * ptab (:,:,jl) 
    96618      END DO 
    97619       
    98620   END FUNCTION fice_cell_ave 
    99621    
    100    FUNCTION fice_ice_ave ( ptab) 
     622    
     623   FUNCTION fice_ice_ave ( ptab ) 
    101624      !!-------------------------------------------------------------------------- 
    102625      !! * Compute average over categories, for ice covered part of grid cell 
     
    106629 
    107630      fice_ice_ave (:,:) = 0.0_wp 
    108       WHERE ( at_i (:,:) .GT. 0.0_wp ) fice_ice_ave (:,:) = fice_cell_ave ( ptab (:,:,:)) / at_i (:,:) 
     631      WHERE ( at_i (:,:) > 0.0_wp ) fice_ice_ave (:,:) = fice_cell_ave ( ptab (:,:,:)) / at_i (:,:) 
    109632 
    110633   END FUNCTION fice_ice_ave 
    111634 
    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 !!gm : remark   old_...   should becomes ...b  as tn versus tb   
    293          old_a_i  (:,:,:)   = a_i  (:,:,:)     ! ice area 
    294          old_e_i  (:,:,:,:) = e_i  (:,:,:,:)   ! ice thermal energy 
    295          old_v_i  (:,:,:)   = v_i  (:,:,:)     ! ice volume 
    296          old_v_s  (:,:,:)   = v_s  (:,:,:)     ! snow volume  
    297          old_e_s  (:,:,:,:) = e_s  (:,:,:,:)   ! snow thermal energy 
    298          old_smv_i(:,:,:)   = smv_i(:,:,:)     ! salt content 
    299          old_oa_i (:,:,:)   = oa_i (:,:,:)     ! areal age content 
    300          old_u_ice(:,:)     = u_ice(:,:) 
    301          old_v_ice(:,:)     = v_ice(:,:) 
    302  
    303          ! trends    !!gm is it truly necessary ??? 
    304          d_a_i_thd  (:,:,:)   = 0._wp   ;   d_a_i_trp  (:,:,:)   = 0._wp 
    305          d_v_i_thd  (:,:,:)   = 0._wp   ;   d_v_i_trp  (:,:,:)   = 0._wp 
    306          d_e_i_thd  (:,:,:,:) = 0._wp   ;   d_e_i_trp  (:,:,:,:) = 0._wp 
    307          d_v_s_thd  (:,:,:)   = 0._wp   ;   d_v_s_trp  (:,:,:)   = 0._wp 
    308          d_e_s_thd  (:,:,:,:) = 0._wp   ;   d_e_s_trp  (:,:,:,:) = 0._wp 
    309          d_smv_i_thd(:,:,:)   = 0._wp   ;   d_smv_i_trp(:,:,:)   = 0._wp 
    310          d_oa_i_thd (:,:,:)   = 0._wp   ;   d_oa_i_trp (:,:,:)   = 0._wp 
    311          d_u_ice_dyn(:,:)     = 0._wp   ;   d_v_ice_dyn(:,:)     = 0._wp 
    312  
    313          ! salt, heat and mass fluxes 
    314          sfx    (:,:) = 0._wp   ; 
    315          sfx_bri(:,:) = 0._wp   ;   sfx_dyn(:,:) = 0._wp  
    316          sfx_sni(:,:) = 0._wp   ;   sfx_opw(:,:) = 0._wp 
    317          sfx_bog(:,:) = 0._wp   ;   sfx_dyn(:,:) = 0._wp 
    318          sfx_bom(:,:) = 0._wp   ;   sfx_sum(:,:) = 0._wp 
    319          sfx_res(:,:) = 0._wp 
    320  
    321          wfx_snw(:,:) = 0._wp   ;   wfx_ice(:,:) = 0._wp 
    322          wfx_sni(:,:) = 0._wp   ;   wfx_opw(:,:) = 0._wp 
    323          wfx_bog(:,:) = 0._wp   ;   wfx_dyn(:,:) = 0._wp 
    324          wfx_bom(:,:) = 0._wp   ;   wfx_sum(:,:) = 0._wp 
    325          wfx_res(:,:) = 0._wp   ;   wfx_sub(:,:) = 0._wp 
    326          wfx_spr(:,:) = 0._wp   ;    
    327  
    328          hfx_in (:,:) = 0._wp   ;   hfx_out(:,:) = 0._wp 
    329          hfx_thd(:,:) = 0._wp   ;    
    330          hfx_snw(:,:) = 0._wp   ;   hfx_opw(:,:) = 0._wp 
    331          hfx_bog(:,:) = 0._wp   ;   hfx_dyn(:,:) = 0._wp 
    332          hfx_bom(:,:) = 0._wp   ;   hfx_sum(:,:) = 0._wp 
    333          hfx_res(:,:) = 0._wp   ;   hfx_sub(:,:) = 0._wp 
    334          hfx_spr(:,:) = 0._wp   ;   hfx_dif(:,:) = 0._wp  
    335          hfx_err(:,:) = 0._wp   ;   hfx_err_rem(:,:) = 0._wp 
    336  
    337          ! 
    338          fhld  (:,:)    = 0._wp  
    339          fmmflx(:,:)    = 0._wp      
    340          ! part of solar radiation transmitted through the ice 
    341          ftr_ice(:,:,:) = 0._wp 
    342  
    343          ! diags 
    344          diag_trp_vi  (:,:) = 0._wp  ; diag_trp_vs(:,:) = 0._wp  ;  diag_trp_ei(:,:) = 0._wp  ;  diag_trp_es(:,:) = 0._wp 
    345          diag_heat_dhc(:,:) = 0._wp   
    346  
    347          ! dynamical invariants 
    348          delta_i(:,:) = 0._wp       ;   divu_i(:,:) = 0._wp       ;   shear_i(:,:) = 0._wp 
    349  
    350                           CALL lim_rst_opn( kt )     ! Open Ice restart file 
    351          ! 
    352          IF( ln_nicep )   CALL lim_prt_state( kt, jiindx, jjindx, 1, ' - Beginning the time step - ' )   ! control print 
    353          ! ---------------------------------------------- 
    354          ! ice dynamics and transport (except in 1D case) 
    355          ! ---------------------------------------------- 
    356          IF( .NOT. lk_c1d ) THEN 
    357                           CALL lim_dyn( kt )              ! Ice dynamics    ( rheology/dynamics ) 
    358                           CALL lim_trp( kt )              ! Ice transport   ( Advection/diffusion ) 
    359                           CALL lim_var_glo2eqv            ! equivalent variables, requested for rafting 
    360          IF( ln_nicep )   CALL lim_prt_state( kt, jiindx, jjindx,-1, ' - ice dyn & trp - ' )   ! control print 
    361                           CALL lim_itd_me                 ! Mechanical redistribution ! (ridging/rafting) 
    362                           CALL lim_var_agg( 1 )  
    363 #if defined key_bdy 
    364                           ! bdy ice thermo  
    365                           CALL lim_var_glo2eqv            ! equivalent variables 
    366                           CALL bdy_ice_lim( kt ) 
    367                           CALL lim_itd_me_zapsmall 
    368                           CALL lim_var_agg(1) 
    369          IF( ln_nicep )   CALL lim_prt_state( kt, jiindx, jjindx, 1, ' - ice thermo bdy - ' )   ! control print 
    370 #endif 
    371                           CALL lim_update1 
    372          ENDIF 
    373 !                         !- Change old values for new values 
    374                           old_u_ice(:,:)   = u_ice (:,:) 
    375                           old_v_ice(:,:)   = v_ice (:,:) 
    376                           old_a_i(:,:,:)   = a_i (:,:,:) 
    377                           old_v_s(:,:,:)   = v_s (:,:,:) 
    378                           old_v_i(:,:,:)   = v_i (:,:,:) 
    379                           old_e_s(:,:,:,:) = e_s (:,:,:,:) 
    380                           old_e_i(:,:,:,:) = e_i (:,:,:,:) 
    381                           old_oa_i(:,:,:)  = oa_i(:,:,:) 
    382                           old_smv_i(:,:,:) = smv_i (:,:,:) 
    383   
    384          ! ---------------------------------------------- 
    385          ! ice thermodynamic 
    386          ! ---------------------------------------------- 
    387                           CALL lim_var_glo2eqv            ! equivalent variables 
    388                           CALL lim_var_agg(1)             ! aggregate ice categories 
    389                           ! previous lead fraction and ice volume for flux calculations 
    390                           pfrld(:,:)   = 1._wp - at_i(:,:) 
    391                           phicif(:,:)  = vt_i(:,:) 
    392                           ! 
    393                           CALL lim_var_bv                 ! bulk brine volume (diag) 
    394                           CALL lim_thd( kt )              ! Ice thermodynamics  
    395                           zcoef = rdt_ice /rday           !  Ice natural aging 
    396                           oa_i(:,:,:) = oa_i(:,:,:) + a_i(:,:,:) * zcoef 
    397          IF( ln_nicep )   CALL lim_prt_state( kt, jiindx, jjindx, 1, ' - ice thermodyn. - ' )   ! control print 
    398                           CALL lim_itd_th( kt )           !  Remap ice categories, lateral accretion  ! 
    399                           CALL lim_var_agg( 1 )           ! requested by limupdate 
    400                           CALL lim_update2                ! Global variables update 
    401  
    402                           CALL lim_var_glo2eqv            ! equivalent variables (outputs) 
    403                           CALL lim_var_agg(2)             ! aggregate ice thickness categories 
    404          IF( ln_nicep )   CALL lim_prt_state( kt, jiindx, jjindx, 2, ' - Final state - ' )   ! control print 
    405          ! 
    406                           CALL lim_sbc_flx( kt )     ! Update surface ocean mass, heat and salt fluxes 
    407          ! 
    408          IF( ln_nicep )   CALL lim_prt_state( kt, jiindx, jjindx, 3, ' - Final state lim_sbc - ' )   ! control print 
    409          ! 
    410          !                                           ! Diagnostics and outputs  
    411          IF (ln_limdiaout) CALL lim_diahsb 
    412  
    413                           CALL lim_wri( 1  )              ! Ice outputs  
    414  
    415          IF( kt == nit000 .AND. ln_rstart )   & 
    416             &             CALL iom_close( numrir )        ! clem: close input ice restart file 
    417          ! 
    418          IF( lrst_ice )   CALL lim_rst_write( kt )        ! Ice restart file  
    419                           CALL lim_var_glo2eqv            ! ??? 
    420          ! 
    421          IF( ln_nicep )   CALL lim_ctl( kt )              ! alerts in case of model crash 
    422          ! 
    423       ENDIF                                    ! End sea-ice time step only 
    424  
    425       !                                        !--------------------------! 
    426       !                                        !  at all ocean time step  ! 
    427       !                                        !--------------------------! 
    428       !                                                
    429       !                                              ! Update surface ocean stresses (only in ice-dynamic case) 
    430       !                                                   ! otherwise the atm.-ocean stresses are used everywhere 
    431       IF( ln_limdyn )     CALL lim_sbc_tau( kt, ub(:,:,1), vb(:,:,1) )  ! using before instantaneous surf. currents 
    432        
    433 !!gm   remark, the ocean-ice stress is not saved in ice diag call above .....  find a solution!!! 
    434       CALL wrk_dealloc( jpi,jpj,jpl, zalb_ice_os, zalb_ice_cs, zalb_ice ) 
    435  
    436       IF( lk_cpl ) THEN 
    437          IF ( ln_iceflx_ave .OR. ln_iceflx_linear ) & 
    438             &    CALL wrk_dealloc( jpi, jpj, ztem_ice_all , zalb_ice_all , z_qsr_ice_all, z_qns_ice_all,   & 
    439             &                                z_qla_ice_all, z_dqns_ice_all, z_dqla_ice_all) 
    440       ENDIF 
    441       ! 
    442       IF( nn_timing == 1 )  CALL timing_stop('sbc_ice_lim') 
    443       ! 
    444    END SUBROUTINE sbc_ice_lim 
    445  
    446  
    447    SUBROUTINE lim_ctl( kt ) 
    448       !!----------------------------------------------------------------------- 
    449       !!                   ***  ROUTINE lim_ctl ***  
    450       !!                  
    451       !! ** Purpose :   Alerts in case of model crash 
    452       !!------------------------------------------------------------------- 
    453       INTEGER, INTENT(in) ::   kt      ! ocean time step 
    454       INTEGER  ::   ji, jj, jk,  jl   ! dummy loop indices 
    455       INTEGER  ::   inb_altests       ! number of alert tests (max 20) 
    456       INTEGER  ::   ialert_id         ! number of the current alert 
    457       REAL(wp) ::   ztmelts           ! ice layer melting point 
    458       CHARACTER (len=30), DIMENSION(20)      ::   cl_alname   ! name of alert 
    459       INTEGER           , DIMENSION(20)      ::   inb_alp     ! number of alerts positive 
    460       !!------------------------------------------------------------------- 
    461  
    462       inb_altests = 10 
    463       inb_alp(:)  =  0 
    464  
    465       ! Alert if incompatible volume and concentration 
    466       ialert_id = 2 ! reference number of this alert 
    467       cl_alname(ialert_id) = ' Incompat vol and con         '    ! name of the alert 
    468  
    469       DO jl = 1, jpl 
    470          DO jj = 1, jpj 
    471             DO ji = 1, jpi 
    472                IF(  v_i(ji,jj,jl) /= 0._wp   .AND.   a_i(ji,jj,jl) == 0._wp   ) THEN 
    473                   !WRITE(numout,*) ' ALERTE 2 :   Incompatible volume and concentration ' 
    474                   !WRITE(numout,*) ' at_i     ', at_i(ji,jj) 
    475                   !WRITE(numout,*) ' Point - category', ji, jj, jl 
    476                   !WRITE(numout,*) ' a_i *** a_i_old ', a_i      (ji,jj,jl), old_a_i  (ji,jj,jl) 
    477                   !WRITE(numout,*) ' v_i *** v_i_old ', v_i      (ji,jj,jl), old_v_i  (ji,jj,jl) 
    478                   !WRITE(numout,*) ' d_a_i_thd/trp   ', d_a_i_thd(ji,jj,jl), d_a_i_trp(ji,jj,jl) 
    479                   !WRITE(numout,*) ' d_v_i_thd/trp   ', d_v_i_thd(ji,jj,jl), d_v_i_trp(ji,jj,jl) 
    480                   inb_alp(ialert_id) = inb_alp(ialert_id) + 1 
    481                ENDIF 
    482             END DO 
    483          END DO 
    484       END DO 
    485  
    486       ! Alerte if very thick ice 
    487       ialert_id = 3 ! reference number of this alert 
    488       cl_alname(ialert_id) = ' Very thick ice               ' ! name of the alert 
    489       jl = jpl  
    490       DO jj = 1, jpj 
    491          DO ji = 1, jpi 
    492             IF(   ht_i(ji,jj,jl)  >  50._wp   ) THEN 
    493                !CALL lim_prt_state( kt, ji, jj, 2, ' ALERTE 3 :   Very thick ice ' ) 
    494                inb_alp(ialert_id) = inb_alp(ialert_id) + 1 
    495             ENDIF 
    496          END DO 
    497       END DO 
    498  
    499       ! Alert if very fast ice 
    500       ialert_id = 4 ! reference number of this alert 
    501       cl_alname(ialert_id) = ' Very fast ice               ' ! name of the alert 
    502       DO jj = 1, jpj 
    503          DO ji = 1, jpi 
    504             IF(   MAX( ABS( u_ice(ji,jj) ), ABS( v_ice(ji,jj) ) ) > 1.5  .AND.  & 
    505                &  at_i(ji,jj) > 0._wp   ) THEN 
    506                !CALL lim_prt_state( kt, ji, jj, 1, ' ALERTE 4 :   Very fast ice ' ) 
    507                !WRITE(numout,*) ' ice strength             : ', strength(ji,jj) 
    508                !WRITE(numout,*) ' oceanic stress utau      : ', utau(ji,jj)  
    509                !WRITE(numout,*) ' oceanic stress vtau      : ', vtau(ji,jj) 
    510                !WRITE(numout,*) ' sea-ice stress utau_ice  : ', utau_ice(ji,jj)  
    511                !WRITE(numout,*) ' sea-ice stress vtau_ice  : ', vtau_ice(ji,jj) 
    512                !WRITE(numout,*) ' oceanic speed u          : ', u_oce(ji,jj) 
    513                !WRITE(numout,*) ' oceanic speed v          : ', v_oce(ji,jj) 
    514                !WRITE(numout,*) ' sst                      : ', sst_m(ji,jj) 
    515                !WRITE(numout,*) ' sss                      : ', sss_m(ji,jj) 
    516                !WRITE(numout,*)  
    517                inb_alp(ialert_id) = inb_alp(ialert_id) + 1 
    518             ENDIF 
    519          END DO 
    520       END DO 
    521  
    522       ! Alert if there is ice on continents 
    523       ialert_id = 6 ! reference number of this alert 
    524       cl_alname(ialert_id) = ' Ice on continents           ' ! name of the alert 
    525       DO jj = 1, jpj 
    526          DO ji = 1, jpi 
    527             IF(   tms(ji,jj) <= 0._wp   .AND.   at_i(ji,jj) > 0._wp   ) THEN  
    528                !CALL lim_prt_state( kt, ji, jj, 1, ' ALERTE 6 :   Ice on continents ' ) 
    529                !WRITE(numout,*) ' masks s, u, v        : ', tms(ji,jj), tmu(ji,jj), tmv(ji,jj)  
    530                !WRITE(numout,*) ' sst                  : ', sst_m(ji,jj) 
    531                !WRITE(numout,*) ' sss                  : ', sss_m(ji,jj) 
    532                !WRITE(numout,*) ' at_i(ji,jj)          : ', at_i(ji,jj) 
    533                !WRITE(numout,*) ' v_ice(ji,jj)         : ', v_ice(ji,jj) 
    534                !WRITE(numout,*) ' v_ice(ji,jj-1)       : ', v_ice(ji,jj-1) 
    535                !WRITE(numout,*) ' u_ice(ji-1,jj)       : ', u_ice(ji-1,jj) 
    536                !WRITE(numout,*) ' u_ice(ji,jj)         : ', v_ice(ji,jj) 
    537                ! 
    538                inb_alp(ialert_id) = inb_alp(ialert_id) + 1 
    539             ENDIF 
    540          END DO 
    541       END DO 
    542  
    543 ! 
    544 !     ! Alert if very fresh ice 
    545       ialert_id = 7 ! reference number of this alert 
    546       cl_alname(ialert_id) = ' Very fresh ice               ' ! name of the alert 
    547       DO jl = 1, jpl 
    548          DO jj = 1, jpj 
    549             DO ji = 1, jpi 
    550                IF( sm_i(ji,jj,jl) < 0.1 .AND. a_i(ji,jj,jl) > 0._wp ) THEN 
    551 !                 CALL lim_prt_state(kt,ji,jj,1, ' ALERTE 7 :   Very fresh ice ' ) 
    552 !                 WRITE(numout,*) ' sst                  : ', sst_m(ji,jj) 
    553 !                 WRITE(numout,*) ' sss                  : ', sss_m(ji,jj) 
    554 !                 WRITE(numout,*)  
    555                   inb_alp(ialert_id) = inb_alp(ialert_id) + 1 
    556                ENDIF 
    557             END DO 
    558          END DO 
    559       END DO 
    560 ! 
    561  
    562 !     ! Alert if too old ice 
    563       ialert_id = 9 ! reference number of this alert 
    564       cl_alname(ialert_id) = ' Very old   ice               ' ! name of the alert 
    565       DO jl = 1, jpl 
    566          DO jj = 1, jpj 
    567             DO ji = 1, jpi 
    568                IF ( ( ( ABS( o_i(ji,jj,jl) ) > rdt_ice ) .OR. & 
    569                       ( ABS( o_i(ji,jj,jl) ) < 0._wp) ) .AND. & 
    570                              ( a_i(ji,jj,jl) > 0._wp ) ) THEN 
    571                   !CALL lim_prt_state( kt, ji, jj, 1, ' ALERTE 9 :   Wrong ice age ') 
    572                   inb_alp(ialert_id) = inb_alp(ialert_id) + 1 
    573                ENDIF 
    574             END DO 
    575          END DO 
    576       END DO 
    577   
    578       ! Alert on salt flux 
    579       ialert_id = 5 ! reference number of this alert 
    580       cl_alname(ialert_id) = ' High salt flux               ' ! name of the alert 
    581       DO jj = 1, jpj 
    582          DO ji = 1, jpi 
    583             IF( ABS( sfx (ji,jj) ) .GT. 1.0e-2 ) THEN  ! = 1 psu/day for 1m ocean depth 
    584                !CALL lim_prt_state( kt, ji, jj, 3, ' ALERTE 5 :   High salt flux ' ) 
    585                !DO jl = 1, jpl 
    586                   !WRITE(numout,*) ' Category no: ', jl 
    587                   !WRITE(numout,*) ' a_i        : ', a_i      (ji,jj,jl) , ' old_a_i    : ', old_a_i  (ji,jj,jl)    
    588                   !WRITE(numout,*) ' d_a_i_trp  : ', d_a_i_trp(ji,jj,jl) , ' d_a_i_thd  : ', d_a_i_thd(ji,jj,jl)  
    589                   !WRITE(numout,*) ' v_i        : ', v_i      (ji,jj,jl) , ' old_v_i    : ', old_v_i  (ji,jj,jl)    
    590                   !WRITE(numout,*) ' d_v_i_trp  : ', d_v_i_trp(ji,jj,jl) , ' d_v_i_thd  : ', d_v_i_thd(ji,jj,jl)  
    591                   !WRITE(numout,*) ' ' 
    592                !END DO 
    593                inb_alp(ialert_id) = inb_alp(ialert_id) + 1 
    594             ENDIF 
    595          END DO 
    596       END DO 
    597  
    598       ! Alert if qns very big 
    599       ialert_id = 8 ! reference number of this alert 
    600       cl_alname(ialert_id) = ' fnsolar very big             ' ! name of the alert 
    601       DO jj = 1, jpj 
    602          DO ji = 1, jpi 
    603             IF( ABS( qns(ji,jj) ) > 1500._wp .AND. at_i(ji,jj) > 0._wp ) THEN 
    604                ! 
    605                !WRITE(numout,*) ' ALERTE 8 :   Very high non-solar heat flux' 
    606                !WRITE(numout,*) ' ji, jj    : ', ji, jj 
    607                !WRITE(numout,*) ' qns       : ', qns(ji,jj) 
    608                !WRITE(numout,*) ' sst       : ', sst_m(ji,jj) 
    609                !WRITE(numout,*) ' sss       : ', sss_m(ji,jj) 
    610                ! 
    611                !CALL lim_prt_state( kt, ji, jj, 2, '   ') 
    612                inb_alp(ialert_id) = inb_alp(ialert_id) + 1 
    613                ! 
    614             ENDIF 
    615          END DO 
    616       END DO 
    617       !+++++ 
    618   
    619       ! Alert if very warm ice 
    620       ialert_id = 10 ! reference number of this alert 
    621       cl_alname(ialert_id) = ' Very warm ice                ' ! name of the alert 
    622       inb_alp(ialert_id) = 0 
    623       DO jl = 1, jpl 
    624          DO jk = 1, nlay_i 
    625             DO jj = 1, jpj 
    626                DO ji = 1, jpi 
    627                   ztmelts    =  -tmut * s_i(ji,jj,jk,jl) + rtt 
    628                   IF( t_i(ji,jj,jk,jl) >= ztmelts  .AND.  v_i(ji,jj,jl) > 1.e-10   & 
    629                      &                             .AND.  a_i(ji,jj,jl) > 0._wp   ) THEN 
    630                      !WRITE(numout,*) ' ALERTE 10 :   Very warm ice' 
    631                      !WRITE(numout,*) ' ji, jj, jk, jl : ', ji, jj, jk, jl 
    632                      !WRITE(numout,*) ' t_i : ', t_i(ji,jj,jk,jl) 
    633                      !WRITE(numout,*) ' e_i : ', e_i(ji,jj,jk,jl) 
    634                      !WRITE(numout,*) ' s_i : ', s_i(ji,jj,jk,jl) 
    635                      !WRITE(numout,*) ' ztmelts : ', ztmelts 
    636                      inb_alp(ialert_id) = inb_alp(ialert_id) + 1 
    637                   ENDIF 
    638                END DO 
    639             END DO 
    640          END DO 
    641       END DO 
    642  
    643       ! sum of the alerts on all processors 
    644       IF( lk_mpp ) THEN 
    645          DO ialert_id = 1, inb_altests 
    646             CALL mpp_sum(inb_alp(ialert_id)) 
    647          END DO 
    648       ENDIF 
    649  
    650       ! print alerts 
    651       IF( lwp ) THEN 
    652          ialert_id = 1                                 ! reference number of this alert 
    653          cl_alname(ialert_id) = ' NO alerte 1      '   ! name of the alert 
    654          WRITE(numout,*) ' time step ',kt 
    655          WRITE(numout,*) ' All alerts at the end of ice model ' 
    656          DO ialert_id = 1, inb_altests 
    657             WRITE(numout,*) ialert_id, cl_alname(ialert_id)//' : ', inb_alp(ialert_id), ' times ! ' 
    658          END DO 
    659       ENDIF 
    660      ! 
    661    END SUBROUTINE lim_ctl 
    662   
    663     
    664    SUBROUTINE lim_prt_state( kt, ki, kj, kn, cd1 ) 
    665       !!----------------------------------------------------------------------- 
    666       !!                   ***  ROUTINE lim_prt_state ***  
    667       !!                  
    668       !! ** Purpose :   Writes global ice state on the (i,j) point  
    669       !!                in ocean.ouput  
    670       !!                3 possibilities exist  
    671       !!                n = 1/-1 -> simple ice state (plus Mechanical Check if -1) 
    672       !!                n = 2    -> exhaustive state 
    673       !!                n = 3    -> ice/ocean salt fluxes 
    674       !! 
    675       !! ** input   :   point coordinates (i,j)  
    676       !!                n : number of the option 
    677       !!------------------------------------------------------------------- 
    678       INTEGER         , INTENT(in) ::   kt      ! ocean time step 
    679       INTEGER         , INTENT(in) ::   ki, kj, kn    ! ocean gridpoint indices 
    680       CHARACTER(len=*), INTENT(in) ::   cd1           ! 
    681       !! 
    682       INTEGER :: jl, ji, jj 
    683       !!------------------------------------------------------------------- 
    684  
    685       DO ji = mi0(ki), mi1(ki) 
    686          DO jj = mj0(kj), mj1(kj) 
    687  
    688             WRITE(numout,*) ' time step ',kt,' ',cd1             ! print title 
    689  
    690             !---------------- 
    691             !  Simple state 
    692             !---------------- 
    693              
    694             IF ( kn == 1 .OR. kn == -1 ) THEN 
    695                WRITE(numout,*) ' lim_prt_state - Point : ',ji,jj 
    696                WRITE(numout,*) ' ~~~~~~~~~~~~~~ ' 
    697                WRITE(numout,*) ' Simple state ' 
    698                WRITE(numout,*) ' masks s,u,v   : ', tms(ji,jj), tmu(ji,jj), tmv(ji,jj) 
    699                WRITE(numout,*) ' lat - long    : ', gphit(ji,jj), glamt(ji,jj) 
    700                WRITE(numout,*) ' Time step     : ', numit 
    701                WRITE(numout,*) ' - Ice drift   ' 
    702                WRITE(numout,*) '   ~~~~~~~~~~~ ' 
    703                WRITE(numout,*) ' u_ice(i-1,j)  : ', u_ice(ji-1,jj) 
    704                WRITE(numout,*) ' u_ice(i  ,j)  : ', u_ice(ji,jj) 
    705                WRITE(numout,*) ' v_ice(i  ,j-1): ', v_ice(ji,jj-1) 
    706                WRITE(numout,*) ' v_ice(i  ,j)  : ', v_ice(ji,jj) 
    707                WRITE(numout,*) ' strength      : ', strength(ji,jj) 
    708                WRITE(numout,*) 
    709                WRITE(numout,*) ' - Cell values ' 
    710                WRITE(numout,*) '   ~~~~~~~~~~~ ' 
    711                WRITE(numout,*) ' cell area     : ', area(ji,jj) 
    712                WRITE(numout,*) ' at_i          : ', at_i(ji,jj)        
    713                WRITE(numout,*) ' vt_i          : ', vt_i(ji,jj)        
    714                WRITE(numout,*) ' vt_s          : ', vt_s(ji,jj)        
    715                DO jl = 1, jpl 
    716                   WRITE(numout,*) ' - Category (', jl,')' 
    717                   WRITE(numout,*) ' a_i           : ', a_i(ji,jj,jl) 
    718                   WRITE(numout,*) ' ht_i          : ', ht_i(ji,jj,jl) 
    719                   WRITE(numout,*) ' ht_s          : ', ht_s(ji,jj,jl) 
    720                   WRITE(numout,*) ' v_i           : ', v_i(ji,jj,jl) 
    721                   WRITE(numout,*) ' v_s           : ', v_s(ji,jj,jl) 
    722                   WRITE(numout,*) ' e_s           : ', e_s(ji,jj,1,jl)/1.0e9 
    723                   WRITE(numout,*) ' e_i           : ', e_i(ji,jj,1:nlay_i,jl)/1.0e9 
    724                   WRITE(numout,*) ' t_su          : ', t_su(ji,jj,jl) 
    725                   WRITE(numout,*) ' t_snow        : ', t_s(ji,jj,1,jl) 
    726                   WRITE(numout,*) ' t_i           : ', t_i(ji,jj,1:nlay_i,jl) 
    727                   WRITE(numout,*) ' sm_i          : ', sm_i(ji,jj,jl) 
    728                   WRITE(numout,*) ' smv_i         : ', smv_i(ji,jj,jl) 
    729                   WRITE(numout,*) 
    730                END DO 
    731             ENDIF 
    732             IF( kn == -1 ) THEN 
    733                WRITE(numout,*) ' Mechanical Check ************** ' 
    734                WRITE(numout,*) ' Check what means ice divergence ' 
    735                WRITE(numout,*) ' Total ice concentration ', at_i (ji,jj) 
    736                WRITE(numout,*) ' Total lead fraction     ', ato_i(ji,jj) 
    737                WRITE(numout,*) ' Sum of both             ', ato_i(ji,jj) + at_i(ji,jj) 
    738                WRITE(numout,*) ' Sum of both minus 1     ', ato_i(ji,jj) + at_i(ji,jj) - 1.00 
    739             ENDIF 
    740              
    741  
    742             !-------------------- 
    743             !  Exhaustive state 
    744             !-------------------- 
    745              
    746             IF ( kn .EQ. 2 ) THEN 
    747                WRITE(numout,*) ' lim_prt_state - Point : ',ji,jj 
    748                WRITE(numout,*) ' ~~~~~~~~~~~~~~ ' 
    749                WRITE(numout,*) ' Exhaustive state ' 
    750                WRITE(numout,*) ' lat - long ', gphit(ji,jj), glamt(ji,jj) 
    751                WRITE(numout,*) ' Time step ', numit 
    752                WRITE(numout,*)  
    753                WRITE(numout,*) ' - Cell values ' 
    754                WRITE(numout,*) '   ~~~~~~~~~~~ ' 
    755                WRITE(numout,*) ' cell area     : ', area(ji,jj) 
    756                WRITE(numout,*) ' at_i          : ', at_i(ji,jj)        
    757                WRITE(numout,*) ' vt_i          : ', vt_i(ji,jj)        
    758                WRITE(numout,*) ' vt_s          : ', vt_s(ji,jj)        
    759                WRITE(numout,*) ' u_ice(i-1,j)  : ', u_ice(ji-1,jj) 
    760                WRITE(numout,*) ' u_ice(i  ,j)  : ', u_ice(ji,jj) 
    761                WRITE(numout,*) ' v_ice(i  ,j-1): ', v_ice(ji,jj-1) 
    762                WRITE(numout,*) ' v_ice(i  ,j)  : ', v_ice(ji,jj) 
    763                WRITE(numout,*) ' strength      : ', strength(ji,jj) 
    764                WRITE(numout,*) ' d_u_ice_dyn   : ', d_u_ice_dyn(ji,jj), ' d_v_ice_dyn   : ', d_v_ice_dyn(ji,jj) 
    765                WRITE(numout,*) ' old_u_ice     : ', old_u_ice(ji,jj)  , ' old_v_ice     : ', old_v_ice(ji,jj)   
    766                WRITE(numout,*) 
    767                 
    768                DO jl = 1, jpl 
    769                   WRITE(numout,*) ' - Category (',jl,')' 
    770                   WRITE(numout,*) '   ~~~~~~~~         '  
    771                   WRITE(numout,*) ' ht_i       : ', ht_i(ji,jj,jl)             , ' ht_s       : ', ht_s(ji,jj,jl) 
    772                   WRITE(numout,*) ' t_i        : ', t_i(ji,jj,1:nlay_i,jl) 
    773                   WRITE(numout,*) ' t_su       : ', t_su(ji,jj,jl)             , ' t_s        : ', t_s(ji,jj,1,jl) 
    774                   WRITE(numout,*) ' sm_i       : ', sm_i(ji,jj,jl)             , ' o_i        : ', o_i(ji,jj,jl) 
    775                   WRITE(numout,*) ' a_i        : ', a_i(ji,jj,jl)              , ' old_a_i    : ', old_a_i(ji,jj,jl)    
    776                   WRITE(numout,*) ' d_a_i_trp  : ', d_a_i_trp(ji,jj,jl)        , ' d_a_i_thd  : ', d_a_i_thd(ji,jj,jl)  
    777                   WRITE(numout,*) ' v_i        : ', v_i(ji,jj,jl)              , ' old_v_i    : ', old_v_i(ji,jj,jl)    
    778                   WRITE(numout,*) ' d_v_i_trp  : ', d_v_i_trp(ji,jj,jl)        , ' d_v_i_thd  : ', d_v_i_thd(ji,jj,jl)  
    779                   WRITE(numout,*) ' v_s        : ', v_s(ji,jj,jl)              , ' old_v_s    : ', old_v_s(ji,jj,jl)   
    780                   WRITE(numout,*) ' d_v_s_trp  : ', d_v_s_trp(ji,jj,jl)        , ' d_v_s_thd  : ', d_v_s_thd(ji,jj,jl) 
    781                   WRITE(numout,*) ' e_i1       : ', e_i(ji,jj,1,jl)/1.0e9      , ' old_ei1    : ', old_e_i(ji,jj,1,jl)/1.0e9  
    782                   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 
    783                   WRITE(numout,*) ' e_i2       : ', e_i(ji,jj,2,jl)/1.0e9      , ' old_ei2    : ', old_e_i(ji,jj,2,jl)/1.0e9   
    784                   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 
    785                   WRITE(numout,*) ' e_snow     : ', e_s(ji,jj,1,jl)            , ' old_e_snow : ', old_e_s(ji,jj,1,jl)  
    786                   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) 
    787                   WRITE(numout,*) ' smv_i      : ', smv_i(ji,jj,jl)            , ' old_smv_i  : ', old_smv_i(ji,jj,jl)    
    788                   WRITE(numout,*) ' d_smv_i_trp: ', d_smv_i_trp(ji,jj,jl)      , ' d_smv_i_thd: ', d_smv_i_thd(ji,jj,jl)  
    789                   WRITE(numout,*) ' oa_i       : ', oa_i(ji,jj,jl)             , ' old_oa_i   : ', old_oa_i(ji,jj,jl) 
    790                   WRITE(numout,*) ' d_oa_i_trp : ', d_oa_i_trp(ji,jj,jl)       , ' d_oa_i_thd : ', d_oa_i_thd(ji,jj,jl) 
    791                END DO !jl 
    792                 
    793                WRITE(numout,*) 
    794                WRITE(numout,*) ' - Heat / FW fluxes ' 
    795                WRITE(numout,*) '   ~~~~~~~~~~~~~~~~ ' 
    796                WRITE(numout,*) ' - Heat fluxes in and out the ice ***' 
    797                WRITE(numout,*) ' qsr_ini       : ', pfrld(ji,jj) * qsr(ji,jj) + SUM( old_a_i(ji,jj,:) * qsr_ice(ji,jj,:) ) 
    798                WRITE(numout,*) ' qns_ini       : ', pfrld(ji,jj) * qns(ji,jj) + SUM( old_a_i(ji,jj,:) * qns_ice(ji,jj,:) ) 
    799                WRITE(numout,*) 
    800                WRITE(numout,*)  
    801                WRITE(numout,*) ' sst        : ', sst_m(ji,jj)   
    802                WRITE(numout,*) ' sss        : ', sss_m(ji,jj)   
    803                WRITE(numout,*)  
    804                WRITE(numout,*) ' - Stresses ' 
    805                WRITE(numout,*) '   ~~~~~~~~ ' 
    806                WRITE(numout,*) ' utau_ice   : ', utau_ice(ji,jj)  
    807                WRITE(numout,*) ' vtau_ice   : ', vtau_ice(ji,jj) 
    808                WRITE(numout,*) ' utau       : ', utau    (ji,jj)  
    809                WRITE(numout,*) ' vtau       : ', vtau    (ji,jj) 
    810                WRITE(numout,*) ' oc. vel. u : ', u_oce   (ji,jj) 
    811                WRITE(numout,*) ' oc. vel. v : ', v_oce   (ji,jj) 
    812             ENDIF 
    813              
    814             !--------------------- 
    815             ! Salt / heat fluxes 
    816             !--------------------- 
    817              
    818             IF ( kn .EQ. 3 ) THEN 
    819                WRITE(numout,*) ' lim_prt_state - Point : ',ji,jj 
    820                WRITE(numout,*) ' ~~~~~~~~~~~~~~ ' 
    821                WRITE(numout,*) ' - Salt / Heat Fluxes ' 
    822                WRITE(numout,*) '   ~~~~~~~~~~~~~~~~ ' 
    823                WRITE(numout,*) ' lat - long ', gphit(ji,jj), glamt(ji,jj) 
    824                WRITE(numout,*) ' Time step ', numit 
    825                WRITE(numout,*) 
    826                WRITE(numout,*) ' - Heat fluxes at bottom interface ***' 
    827                WRITE(numout,*) ' qsr       : ', qsr(ji,jj) 
    828                WRITE(numout,*) ' qns       : ', qns(ji,jj) 
    829                WRITE(numout,*) 
    830                WRITE(numout,*) ' hfx_mass     : ', hfx_thd(ji,jj) + hfx_dyn(ji,jj) + hfx_snw(ji,jj) + hfx_res(ji,jj) 
    831                WRITE(numout,*) ' hfx_in       : ', hfx_in(ji,jj) 
    832                WRITE(numout,*) ' hfx_out      : ', hfx_out(ji,jj) 
    833                WRITE(numout,*) ' dhc          : ', diag_heat_dhc(ji,jj)               
    834                WRITE(numout,*) 
    835                WRITE(numout,*) ' hfx_dyn      : ', hfx_dyn(ji,jj) 
    836                WRITE(numout,*) ' hfx_thd      : ', hfx_thd(ji,jj) 
    837                WRITE(numout,*) ' hfx_res      : ', hfx_res(ji,jj) 
    838                WRITE(numout,*) ' fhtur        : ', fhtur(ji,jj)  
    839                WRITE(numout,*) ' qlead        : ', qlead(ji,jj) * r1_rdtice 
    840                WRITE(numout,*) 
    841                WRITE(numout,*) ' - Salt fluxes at bottom interface ***' 
    842                WRITE(numout,*) ' emp       : ', emp    (ji,jj) 
    843                WRITE(numout,*) ' sfx       : ', sfx    (ji,jj) 
    844                WRITE(numout,*) ' sfx_res   : ', sfx_res(ji,jj) 
    845                WRITE(numout,*) ' sfx_bri   : ', sfx_bri(ji,jj) 
    846                WRITE(numout,*) ' sfx_dyn   : ', sfx_dyn(ji,jj) 
    847                WRITE(numout,*) 
    848                WRITE(numout,*) ' - Momentum fluxes ' 
    849                WRITE(numout,*) ' utau      : ', utau(ji,jj)  
    850                WRITE(numout,*) ' vtau      : ', vtau(ji,jj) 
    851             ENDIF  
    852             WRITE(numout,*) ' ' 
    853             ! 
    854          END DO 
    855       END DO 
    856  
    857    END SUBROUTINE lim_prt_state 
    858635 
    859636#else 
     
    865642      WRITE(*,*) 'sbc_ice_lim: You should not have seen this print! error?', kt, kblk 
    866643   END SUBROUTINE sbc_ice_lim 
     644   SUBROUTINE sbc_lim_init                 ! Dummy routine 
     645   END SUBROUTINE sbc_lim_init 
    867646#endif 
    868647 
Note: See TracChangeset for help on using the changeset viewer.