MODULE sbcice_lim !!====================================================================== !! *** MODULE sbcice_lim *** !! Surface module : update the ocean surface boundary condition over ice !! & covered area using LIM sea-ice model !! Sea-Ice model : LIM-3 Sea ice model time-stepping !!===================================================================== !! History : 2.0 ! 2006-12 (M. Vancoppenolle) Original code !! 3.0 ! 2008-02 (C. Talandier) Surface module from icestp.F90 !! - ! 2008-04 (G. Madec) sltyle and lim_ctl routine !! 3.3 ! 2010-11 (G. Madec) ice-ocean stress always computed at each ocean time-step !! 3.4 ! 2011-01 (A Porter) dynamical allocation !! - ! 2012-10 (C. Rousset) add lim_diahsb !! 3.6 ! 2014-07 (M. Vancoppenolle, G. Madec, O. Marti) revise coupled interface !!---------------------------------------------------------------------- #if defined key_lim3 !!---------------------------------------------------------------------- !! 'key_lim3' : LIM 3.0 sea-ice model !!---------------------------------------------------------------------- !! sbc_ice_lim : sea-ice model time-stepping and update ocean sbc over ice-covered area !! lim_ctl : alerts in case of ice model crash !! lim_prt_state : ice control print at a given grid point !!---------------------------------------------------------------------- USE oce ! ocean dynamics and tracers USE dom_oce ! ocean space and time domain USE par_ice ! sea-ice parameters USE ice ! LIM-3: ice variables USE iceini ! LIM-3: ice initialisation USE dom_ice ! LIM-3: ice domain USE sbc_oce ! Surface boundary condition: ocean fields USE sbc_ice ! Surface boundary condition: ice fields USE sbcblk_core ! Surface boundary condition: CORE bulk USE sbcblk_clio ! Surface boundary condition: CLIO bulk USE sbccpl ! Surface boundary condition: coupled interface USE albedo ! ocean & ice albedo USE phycst ! Define parameters for the routines USE eosbn2 ! equation of state USE limdyn ! Ice dynamics USE limtrp ! Ice transport USE limthd ! Ice thermodynamics USE limitd_th ! Thermodynamics on ice thickness distribution USE limitd_me ! Mechanics on ice thickness distribution USE limsbc ! sea surface boundary condition USE limdiahsb ! Ice budget diagnostics USE limwri ! Ice outputs USE limrst ! Ice restarts USE limupdate1 ! update of global variables USE limupdate2 ! update of global variables USE limvar ! Ice variables switch USE c1d ! 1D vertical configuration USE lbclnk ! lateral boundary condition - MPP link USE lib_mpp ! MPP library USE wrk_nemo ! work arrays USE timing ! Timing USE iom ! I/O manager library USE in_out_manager ! I/O manager USE prtctl ! Print control USE lib_fortran ! #if defined key_bdy USE bdyice_lim ! unstructured open boundary data (bdy_ice_lim routine) #endif IMPLICIT NONE PRIVATE PUBLIC sbc_ice_lim ! routine called by sbcmod.F90 PUBLIC lim_prt_state !! * Substitutions # include "domzgr_substitute.h90" # include "vectopt_loop_substitute.h90" !!---------------------------------------------------------------------- !! NEMO/OPA 4.0 , UCL NEMO Consortium (2011) !! $Id$ !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) !!---------------------------------------------------------------------- CONTAINS !!====================================================================== SUBROUTINE sbc_ice_lim( kt, kblk ) !!--------------------------------------------------------------------- !! *** ROUTINE sbc_ice_lim *** !! !! ** Purpose : update the ocean surface boundary condition via the !! Louvain la Neuve Sea Ice Model time stepping !! !! ** Method : ice model time stepping !! - call the ice dynamics routine !! - call the ice advection/diffusion routine !! - call the ice thermodynamics routine !! - call the routine that computes mass and !! heat fluxes at the ice/ocean interface !! - save the outputs !! - save the outputs for restart when necessary !! !! ** Action : - time evolution of the LIM sea-ice model !! - update all sbc variables below sea-ice: !! utau, vtau, taum, wndm, qns , qsr, emp , sfx !!--------------------------------------------------------------------- INTEGER, INTENT(in) :: kt ! ocean time step INTEGER, INTENT(in) :: kblk ! type of bulk (=3 CLIO, =4 CORE, =5 COUPLED) !! INTEGER :: jl ! dummy loop index REAL(wp) :: zcoef ! local scalar REAL(wp), POINTER, DIMENSION(:,:,:) :: zalb_os, zalb_cs ! ice albedo under overcast/clear sky REAL(wp), POINTER, DIMENSION(:,:,:) :: zalb_ice ! mean ice albedo (for coupled) !!---------------------------------------------------------------------- IF( nn_timing == 1 ) CALL timing_start('sbc_ice_lim') IF( kt == nit000 ) THEN IF(lwp) WRITE(numout,*) IF(lwp) WRITE(numout,*) 'sbc_ice_lim : update ocean surface boudary condition' IF(lwp) WRITE(numout,*) '~~~~~~~~~~~ via Louvain la Neuve Ice Model (LIM-3) time stepping' ! CALL ice_init ! IF( ln_nicep ) THEN ! control print at a given point jiindx = 15 ; jjindx = 44 IF(lwp) WRITE(numout,*) ' The debugging point is : jiindx : ',jiindx, ' jjindx : ',jjindx ENDIF ENDIF ! !----------------------! IF( MOD( kt-1, nn_fsbc ) == 0 ) THEN ! Ice time-step only ! ! !----------------------! ! ! Bulk Formulae ! ! !----------------! ! u_oce(:,:) = ssu_m(:,:) * umask(:,:,1) ! mean surface ocean current at ice velocity point v_oce(:,:) = ssv_m(:,:) * vmask(:,:,1) ! (C-grid dynamics : U- & V-points as the ocean) ! t_bo(:,:) = ( eos_fzp( sss_m ) + rt0 ) * tmask(:,:,1) + rt0 * ( 1. - tmask(:,:,1) ) ! masked sea surface freezing temperature [Kelvin] ! ! (set to rt0 over land) ! ! Ice albedo CALL wrk_alloc( jpi,jpj,jpl, zalb_os, zalb_cs, zalb_ice ) CALL albedo_ice( t_su, ht_i, ht_s, zalb_cs, zalb_os ) ! cloud-sky and overcast-sky ice albedos SELECT CASE( kblk ) CASE( jp_core , jp_cpl ) ! CORE and COUPLED bulk formulations ! albedo depends on cloud fraction because of non-linear spectral effects zalb_ice(:,:,:) = ( 1. - cldf_ice ) * zalb_cs(:,:,:) + cldf_ice * zalb_os(:,:,:) ! In CLIO the cloud fraction is read in the climatology and the all-sky albedo ! (zalb_ice) is computed within the bulk routine END SELECT ! ! Mask sea ice surface temperature DO jl = 1, jpl t_su(:,:,jl) = t_su(:,:,jl) + rt0 * ( 1. - tmask(:,:,1) ) END DO ! Bulk formulae - provides the following fields: ! utau_ice, vtau_ice : surface ice stress (U- & V-points) [N/m2] ! qsr_ice , qns_ice : solar & non solar heat flux over ice (T-point) [W/m2] ! qla_ice : latent heat flux over ice (T-point) [W/m2] ! dqns_ice, dqla_ice : non solar & latent heat sensistivity (T-point) [W/m2] ! tprecip , sprecip : total & solid precipitation (T-point) [Kg/m2/s] ! fr1_i0 , fr2_i0 : 1sr & 2nd fraction of qsr penetration in ice [%] ! SELECT CASE( kblk ) CASE( jp_clio ) ! CLIO bulk formulation CALL blk_ice_clio( t_su , zalb_cs , zalb_os , zalb_ice , & & utau_ice , vtau_ice , qns_ice , qsr_ice , & & qla_ice , dqns_ice , dqla_ice , & & tprecip , sprecip , & & fr1_i0 , fr2_i0 , cp_ice_msh, jpl ) ! IF( nn_limflx /= 2 ) CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice , & & dqns_ice, qla_ice, dqla_ice, nn_limflx ) CASE( jp_core ) ! CORE bulk formulation CALL blk_ice_core( t_su , u_ice , v_ice , zalb_ice , & & utau_ice , vtau_ice , qns_ice , qsr_ice , & & qla_ice , dqns_ice , dqla_ice , & & tprecip , sprecip , & & fr1_i0 , fr2_i0 , cp_ice_msh, jpl ) ! IF( nn_limflx /= 2 ) CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice , & & dqns_ice, qla_ice, dqla_ice, nn_limflx ) ! CASE ( jp_cpl ) CALL sbc_cpl_ice_tau( utau_ice , vtau_ice ) ! MV -> seb ! CALL sbc_cpl_ice_flx( p_frld=ato_i, palbi=zalb_ice, psst=sst_m, pist=t_su ) ! IF( nn_limflx == 2 ) CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice , & ! & dqns_ice, qla_ice, dqla_ice, nn_limflx ) ! ! Latent heat flux is forced to 0 in coupled : ! ! it is included in qns (non-solar heat flux) ! qla_ice (:,:,:) = 0._wp ! dqla_ice (:,:,:) = 0._wp ! END MV -> seb ! END SELECT ! !----------------------! ! ! LIM-3 time-stepping ! ! !----------------------! ! numit = numit + nn_fsbc ! Ice model time step ! ! ! Store previous ice values a_i_b (:,:,:) = a_i (:,:,:) ! ice area e_i_b (:,:,:,:) = e_i (:,:,:,:) ! ice thermal energy v_i_b (:,:,:) = v_i (:,:,:) ! ice volume v_s_b (:,:,:) = v_s (:,:,:) ! snow volume e_s_b (:,:,:,:) = e_s (:,:,:,:) ! snow thermal energy smv_i_b(:,:,:) = smv_i(:,:,:) ! salt content oa_i_b (:,:,:) = oa_i (:,:,:) ! areal age content u_ice_b(:,:) = u_ice(:,:) v_ice_b(:,:) = v_ice(:,:) ! salt, heat and mass fluxes sfx (:,:) = 0._wp ; sfx_bri(:,:) = 0._wp ; sfx_sni(:,:) = 0._wp ; sfx_opw(:,:) = 0._wp sfx_bog(:,:) = 0._wp ; sfx_dyn(:,:) = 0._wp sfx_bom(:,:) = 0._wp ; sfx_sum(:,:) = 0._wp sfx_res(:,:) = 0._wp wfx_snw(:,:) = 0._wp ; wfx_ice(:,:) = 0._wp wfx_sni(:,:) = 0._wp ; wfx_opw(:,:) = 0._wp wfx_bog(:,:) = 0._wp ; wfx_dyn(:,:) = 0._wp wfx_bom(:,:) = 0._wp ; wfx_sum(:,:) = 0._wp wfx_res(:,:) = 0._wp ; wfx_sub(:,:) = 0._wp wfx_spr(:,:) = 0._wp ; hfx_in (:,:) = 0._wp ; hfx_out(:,:) = 0._wp hfx_thd(:,:) = 0._wp ; hfx_snw(:,:) = 0._wp ; hfx_opw(:,:) = 0._wp hfx_bog(:,:) = 0._wp ; hfx_dyn(:,:) = 0._wp hfx_bom(:,:) = 0._wp ; hfx_sum(:,:) = 0._wp hfx_res(:,:) = 0._wp ; hfx_sub(:,:) = 0._wp hfx_spr(:,:) = 0._wp ; hfx_dif(:,:) = 0._wp hfx_err(:,:) = 0._wp ; hfx_err_rem(:,:) = 0._wp CALL lim_rst_opn( kt ) ! Open Ice restart file ! IF( ln_nicep ) CALL lim_prt_state( kt, jiindx, jjindx, 1, ' - Beginning the time step - ' ) ! control print ! ---------------------------------------------- ! ice dynamics and transport (except in 1D case) ! ---------------------------------------------- IF( .NOT. lk_c1d ) THEN CALL lim_dyn( kt ) ! Ice dynamics ( rheology/dynamics ) CALL lim_trp( kt ) ! Ice transport ( Advection/diffusion ) CALL lim_var_glo2eqv ! equivalent variables, requested for rafting IF( ln_nicep ) CALL lim_prt_state( kt, jiindx, jjindx,-1, ' - ice dyn & trp - ' ) ! control print CALL lim_itd_me ! Mechanical redistribution ! (ridging/rafting) CALL lim_var_agg( 1 ) #if defined key_bdy ! bdy ice thermo CALL lim_var_glo2eqv ! equivalent variables CALL bdy_ice_lim( kt ) CALL lim_itd_me_zapsmall CALL lim_var_agg(1) IF( ln_nicep ) CALL lim_prt_state( kt, jiindx, jjindx, 1, ' - ice thermo bdy - ' ) ! control print #endif CALL lim_update1 ENDIF ! !- Change old values for new values u_ice_b(:,:) = u_ice(:,:) v_ice_b(:,:) = v_ice(:,:) a_i_b (:,:,:) = a_i (:,:,:) v_s_b (:,:,:) = v_s (:,:,:) v_i_b (:,:,:) = v_i (:,:,:) e_s_b (:,:,:,:) = e_s (:,:,:,:) e_i_b (:,:,:,:) = e_i (:,:,:,:) oa_i_b (:,:,:) = oa_i (:,:,:) smv_i_b(:,:,:) = smv_i(:,:,:) ! ---------------------------------------------- ! ice thermodynamic ! ---------------------------------------------- CALL lim_var_glo2eqv ! equivalent variables CALL lim_var_agg(1) ! aggregate ice categories ! previous lead fraction and ice volume for flux calculations pfrld(:,:) = 1._wp - at_i(:,:) phicif(:,:) = vt_i(:,:) ! MV -> seb SELECT CASE( kblk ) CASE ( jp_cpl ) CALL sbc_cpl_ice_flx( p_frld=pfrld, palbi=zalb_ice, psst=sst_m, pist=t_su ) IF( nn_limflx == 2 ) CALL ice_lim_flx( t_su, zalb_ice, qns_ice, qsr_ice , & & dqns_ice, qla_ice, dqla_ice, nn_limflx ) ! Latent heat flux is forced to 0 in coupled : ! it is included in qns (non-solar heat flux) qla_ice (:,:,:) = 0._wp dqla_ice (:,:,:) = 0._wp END SELECT ! END MV -> seb ! CALL lim_var_bv ! bulk brine volume (diag) CALL lim_thd( kt ) ! Ice thermodynamics zcoef = rdt_ice /rday ! Ice natural aging oa_i(:,:,:) = oa_i(:,:,:) + a_i(:,:,:) * zcoef IF( ln_nicep ) CALL lim_prt_state( kt, jiindx, jjindx, 1, ' - ice thermodyn. - ' ) ! control print CALL lim_itd_th( kt ) ! Remap ice categories, lateral accretion ! CALL lim_var_agg( 1 ) ! requested by limupdate CALL lim_update2 ! Global variables update CALL lim_var_glo2eqv ! equivalent variables (outputs) CALL lim_var_agg(2) ! aggregate ice thickness categories IF( ln_nicep ) CALL lim_prt_state( kt, jiindx, jjindx, 2, ' - Final state - ' ) ! control print ! CALL lim_sbc_flx( kt ) ! Update surface ocean mass, heat and salt fluxes ! IF( ln_nicep ) CALL lim_prt_state( kt, jiindx, jjindx, 3, ' - Final state lim_sbc - ' ) ! control print ! ! ! Diagnostics and outputs IF (ln_limdiaout) CALL lim_diahsb CALL lim_wri( 1 ) ! Ice outputs IF( kt == nit000 .AND. ln_rstart ) & & CALL iom_close( numrir ) ! clem: close input ice restart file ! IF( lrst_ice ) CALL lim_rst_write( kt ) ! Ice restart file CALL lim_var_glo2eqv ! ??? ! IF( ln_nicep ) CALL lim_ctl( kt ) ! alerts in case of model crash ! CALL wrk_dealloc( jpi,jpj,jpl, zalb_os, zalb_cs, zalb_ice ) ! ENDIF ! End sea-ice time step only ! !--------------------------! ! ! at all ocean time step ! ! !--------------------------! ! ! ! Update surface ocean stresses (only in ice-dynamic case) ! ! otherwise the atm.-ocean stresses are used everywhere IF( ln_limdyn ) CALL lim_sbc_tau( kt, ub(:,:,1), vb(:,:,1) ) ! using before instantaneous surf. currents !!gm remark, the ocean-ice stress is not saved in ice diag call above ..... find a solution!!! ! IF( nn_timing == 1 ) CALL timing_stop('sbc_ice_lim') ! END SUBROUTINE sbc_ice_lim SUBROUTINE ice_lim_flx( ptn_ice, palb_ice, pqns_ice, pqsr_ice, & & pdqn_ice, pqla_ice, pdql_ice, k_limflx ) !!--------------------------------------------------------------------- !! *** ROUTINE sbc_ice_lim *** !! !! ** Purpose : update the ice surface boundary condition by averaging and / or !! redistributing fluxes on ice categories !! !! ** Method : average then redistribute !! !! ** Action : !!--------------------------------------------------------------------- INTEGER , INTENT(in ) :: k_limflx ! =-1 do nothing; =0 average ; ! =1 average and redistribute ; =2 redistribute REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: ptn_ice ! ice surface temperature REAL(wp), DIMENSION(:,:,:), INTENT(in ) :: palb_ice ! ice albedo REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pqns_ice ! non solar flux REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pqsr_ice ! net solar flux REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pdqn_ice ! non solar flux sensitivity REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pqla_ice ! latent heat flux REAL(wp), DIMENSION(:,:,:), INTENT(inout) :: pdql_ice ! latent heat flux sensitivity ! INTEGER :: jl ! dummy loop index ! REAL(wp), POINTER, DIMENSION(:,:) :: zalb_m ! Mean albedo over all categories REAL(wp), POINTER, DIMENSION(:,:) :: ztem_m ! Mean temperature over all categories ! REAL(wp), POINTER, DIMENSION(:,:) :: z_qsr_m ! Mean solar heat flux over all categories REAL(wp), POINTER, DIMENSION(:,:) :: z_qns_m ! Mean non solar heat flux over all categories REAL(wp), POINTER, DIMENSION(:,:) :: z_qla_m ! Mean latent heat flux over all categories REAL(wp), POINTER, DIMENSION(:,:) :: z_dqn_m ! Mean d(qns)/dT over all categories REAL(wp), POINTER, DIMENSION(:,:) :: z_dql_m ! Mean d(qla)/dT over all categories !!---------------------------------------------------------------------- IF( nn_timing == 1 ) CALL timing_start('ice_lim_flx') ! ! SELECT CASE( k_limflx ) !== averaged on all ice categories ==! CASE( 0 , 1 ) CALL wrk_alloc( jpi,jpj, z_qsr_m, z_qns_m, z_qla_m, z_dqn_m, z_dql_m) ! z_qns_m(:,:) = fice_ice_ave ( pqns_ice (:,:,:) ) z_qsr_m(:,:) = fice_ice_ave ( pqsr_ice (:,:,:) ) z_dqn_m(:,:) = fice_ice_ave ( pdqn_ice (:,:,:) ) z_qla_m(:,:) = fice_ice_ave ( pqla_ice (:,:,:) ) z_dql_m(:,:) = fice_ice_ave ( pdql_ice (:,:,:) ) DO jl = 1, jpl pdqn_ice(:,:,jl) = z_dqn_m(:,:) pdql_ice(:,:,jl) = z_dql_m(:,:) END DO ! DO jl = 1, jpl pqns_ice(:,:,jl) = z_qns_m(:,:) pqsr_ice(:,:,jl) = z_qsr_m(:,:) pqla_ice(:,:,jl) = z_qla_m(:,:) END DO ! CALL wrk_dealloc( jpi,jpj, z_qsr_m, z_qns_m, z_qla_m, z_dqn_m, z_dql_m) END SELECT SELECT CASE( k_limflx ) !== redistribution on all ice categories ==! CASE( 1 , 2 ) CALL wrk_alloc( jpi,jpj, zalb_m, ztem_m ) ! zalb_m(:,:) = fice_ice_ave ( palb_ice (:,:,:) ) ztem_m(:,:) = fice_ice_ave ( ptn_ice (:,:,:) ) DO jl = 1, jpl pqns_ice(:,:,jl) = pqns_ice(:,:,jl) + pdqn_ice(:,:,jl) * (ptn_ice(:,:,jl) - ztem_m(:,:)) pqla_ice(:,:,jl) = pqla_ice(:,:,jl) + pdql_ice(:,:,jl) * (ptn_ice(:,:,jl) - ztem_m(:,:)) pqsr_ice(:,:,jl) = pqsr_ice(:,:,jl) * ( 1._wp - palb_ice(:,:,jl) ) / ( 1._wp - zalb_m(:,:) ) END DO ! CALL wrk_dealloc( jpi,jpj, zalb_m, ztem_m ) END SELECT ! IF( nn_timing == 1 ) CALL timing_stop('ice_lim_flx') ! END SUBROUTINE ice_lim_flx SUBROUTINE lim_ctl( kt ) !!----------------------------------------------------------------------- !! *** ROUTINE lim_ctl *** !! !! ** Purpose : Alerts in case of model crash !!------------------------------------------------------------------- INTEGER, INTENT(in) :: kt ! ocean time step INTEGER :: ji, jj, jk, jl ! dummy loop indices INTEGER :: inb_altests ! number of alert tests (max 20) INTEGER :: ialert_id ! number of the current alert REAL(wp) :: ztmelts ! ice layer melting point CHARACTER (len=30), DIMENSION(20) :: cl_alname ! name of alert INTEGER , DIMENSION(20) :: inb_alp ! number of alerts positive !!------------------------------------------------------------------- inb_altests = 10 inb_alp(:) = 0 ! Alert if incompatible volume and concentration ialert_id = 2 ! reference number of this alert cl_alname(ialert_id) = ' Incompat vol and con ' ! name of the alert DO jl = 1, jpl DO jj = 1, jpj DO ji = 1, jpi IF( v_i(ji,jj,jl) /= 0._wp .AND. a_i(ji,jj,jl) == 0._wp ) THEN !WRITE(numout,*) ' ALERTE 2 : Incompatible volume and concentration ' !WRITE(numout,*) ' at_i ', at_i(ji,jj) !WRITE(numout,*) ' Point - category', ji, jj, jl !WRITE(numout,*) ' a_i *** a_i_b ', a_i (ji,jj,jl), a_i_b (ji,jj,jl) !WRITE(numout,*) ' v_i *** v_i_b ', v_i (ji,jj,jl), v_i_b (ji,jj,jl) !WRITE(numout,*) ' d_a_i_thd/trp ', d_a_i_thd(ji,jj,jl), d_a_i_trp(ji,jj,jl) !WRITE(numout,*) ' d_v_i_thd/trp ', d_v_i_thd(ji,jj,jl), d_v_i_trp(ji,jj,jl) inb_alp(ialert_id) = inb_alp(ialert_id) + 1 ENDIF END DO END DO END DO ! Alerte if very thick ice ialert_id = 3 ! reference number of this alert cl_alname(ialert_id) = ' Very thick ice ' ! name of the alert jl = jpl DO jj = 1, jpj DO ji = 1, jpi IF( ht_i(ji,jj,jl) > 50._wp ) THEN !CALL lim_prt_state( kt, ji, jj, 2, ' ALERTE 3 : Very thick ice ' ) inb_alp(ialert_id) = inb_alp(ialert_id) + 1 ENDIF END DO END DO ! Alert if very fast ice ialert_id = 4 ! reference number of this alert cl_alname(ialert_id) = ' Very fast ice ' ! name of the alert DO jj = 1, jpj DO ji = 1, jpi IF( MAX( ABS( u_ice(ji,jj) ), ABS( v_ice(ji,jj) ) ) > 1.5 .AND. & & at_i(ji,jj) > 0._wp ) THEN !CALL lim_prt_state( kt, ji, jj, 1, ' ALERTE 4 : Very fast ice ' ) !WRITE(numout,*) ' ice strength : ', strength(ji,jj) !WRITE(numout,*) ' oceanic stress utau : ', utau(ji,jj) !WRITE(numout,*) ' oceanic stress vtau : ', vtau(ji,jj) !WRITE(numout,*) ' sea-ice stress utau_ice : ', utau_ice(ji,jj) !WRITE(numout,*) ' sea-ice stress vtau_ice : ', vtau_ice(ji,jj) !WRITE(numout,*) ' oceanic speed u : ', u_oce(ji,jj) !WRITE(numout,*) ' oceanic speed v : ', v_oce(ji,jj) !WRITE(numout,*) ' sst : ', sst_m(ji,jj) !WRITE(numout,*) ' sss : ', sss_m(ji,jj) !WRITE(numout,*) inb_alp(ialert_id) = inb_alp(ialert_id) + 1 ENDIF END DO END DO ! Alert if there is ice on continents ialert_id = 6 ! reference number of this alert cl_alname(ialert_id) = ' Ice on continents ' ! name of the alert DO jj = 1, jpj DO ji = 1, jpi IF( tms(ji,jj) <= 0._wp .AND. at_i(ji,jj) > 0._wp ) THEN !CALL lim_prt_state( kt, ji, jj, 1, ' ALERTE 6 : Ice on continents ' ) !WRITE(numout,*) ' masks s, u, v : ', tms(ji,jj), tmu(ji,jj), tmv(ji,jj) !WRITE(numout,*) ' sst : ', sst_m(ji,jj) !WRITE(numout,*) ' sss : ', sss_m(ji,jj) !WRITE(numout,*) ' at_i(ji,jj) : ', at_i(ji,jj) !WRITE(numout,*) ' v_ice(ji,jj) : ', v_ice(ji,jj) !WRITE(numout,*) ' v_ice(ji,jj-1) : ', v_ice(ji,jj-1) !WRITE(numout,*) ' u_ice(ji-1,jj) : ', u_ice(ji-1,jj) !WRITE(numout,*) ' u_ice(ji,jj) : ', v_ice(ji,jj) ! inb_alp(ialert_id) = inb_alp(ialert_id) + 1 ENDIF END DO END DO ! ! ! Alert if very fresh ice ialert_id = 7 ! reference number of this alert cl_alname(ialert_id) = ' Very fresh ice ' ! name of the alert DO jl = 1, jpl DO jj = 1, jpj DO ji = 1, jpi IF( sm_i(ji,jj,jl) < 0.1 .AND. a_i(ji,jj,jl) > 0._wp ) THEN ! CALL lim_prt_state(kt,ji,jj,1, ' ALERTE 7 : Very fresh ice ' ) ! WRITE(numout,*) ' sst : ', sst_m(ji,jj) ! WRITE(numout,*) ' sss : ', sss_m(ji,jj) ! WRITE(numout,*) inb_alp(ialert_id) = inb_alp(ialert_id) + 1 ENDIF END DO END DO END DO ! ! ! Alert if too old ice ialert_id = 9 ! reference number of this alert cl_alname(ialert_id) = ' Very old ice ' ! name of the alert DO jl = 1, jpl DO jj = 1, jpj DO ji = 1, jpi IF ( ( ( ABS( o_i(ji,jj,jl) ) > rdt_ice ) .OR. & ( ABS( o_i(ji,jj,jl) ) < 0._wp) ) .AND. & ( a_i(ji,jj,jl) > 0._wp ) ) THEN !CALL lim_prt_state( kt, ji, jj, 1, ' ALERTE 9 : Wrong ice age ') inb_alp(ialert_id) = inb_alp(ialert_id) + 1 ENDIF END DO END DO END DO ! Alert on salt flux ialert_id = 5 ! reference number of this alert cl_alname(ialert_id) = ' High salt flux ' ! name of the alert DO jj = 1, jpj DO ji = 1, jpi IF( ABS( sfx (ji,jj) ) .GT. 1.0e-2 ) THEN ! = 1 psu/day for 1m ocean depth !CALL lim_prt_state( kt, ji, jj, 3, ' ALERTE 5 : High salt flux ' ) !DO jl = 1, jpl !WRITE(numout,*) ' Category no: ', jl !WRITE(numout,*) ' a_i : ', a_i (ji,jj,jl) , ' a_i_b : ', a_i_b (ji,jj,jl) !WRITE(numout,*) ' d_a_i_trp : ', d_a_i_trp(ji,jj,jl) , ' d_a_i_thd : ', d_a_i_thd(ji,jj,jl) !WRITE(numout,*) ' v_i : ', v_i (ji,jj,jl) , ' v_i_b : ', v_i_b (ji,jj,jl) !WRITE(numout,*) ' d_v_i_trp : ', d_v_i_trp(ji,jj,jl) , ' d_v_i_thd : ', d_v_i_thd(ji,jj,jl) !WRITE(numout,*) ' ' !END DO inb_alp(ialert_id) = inb_alp(ialert_id) + 1 ENDIF END DO END DO ! Alert if qns very big ialert_id = 8 ! reference number of this alert cl_alname(ialert_id) = ' fnsolar very big ' ! name of the alert DO jj = 1, jpj DO ji = 1, jpi IF( ABS( qns(ji,jj) ) > 1500._wp .AND. at_i(ji,jj) > 0._wp ) THEN ! !WRITE(numout,*) ' ALERTE 8 : Very high non-solar heat flux' !WRITE(numout,*) ' ji, jj : ', ji, jj !WRITE(numout,*) ' qns : ', qns(ji,jj) !WRITE(numout,*) ' sst : ', sst_m(ji,jj) !WRITE(numout,*) ' sss : ', sss_m(ji,jj) ! !CALL lim_prt_state( kt, ji, jj, 2, ' ') inb_alp(ialert_id) = inb_alp(ialert_id) + 1 ! ENDIF END DO END DO !+++++ ! Alert if very warm ice ialert_id = 10 ! reference number of this alert cl_alname(ialert_id) = ' Very warm ice ' ! name of the alert inb_alp(ialert_id) = 0 DO jl = 1, jpl DO jk = 1, nlay_i DO jj = 1, jpj DO ji = 1, jpi ztmelts = -tmut * s_i(ji,jj,jk,jl) + rtt IF( t_i(ji,jj,jk,jl) >= ztmelts .AND. v_i(ji,jj,jl) > 1.e-10 & & .AND. a_i(ji,jj,jl) > 0._wp ) THEN !WRITE(numout,*) ' ALERTE 10 : Very warm ice' !WRITE(numout,*) ' ji, jj, jk, jl : ', ji, jj, jk, jl !WRITE(numout,*) ' t_i : ', t_i(ji,jj,jk,jl) !WRITE(numout,*) ' e_i : ', e_i(ji,jj,jk,jl) !WRITE(numout,*) ' s_i : ', s_i(ji,jj,jk,jl) !WRITE(numout,*) ' ztmelts : ', ztmelts inb_alp(ialert_id) = inb_alp(ialert_id) + 1 ENDIF END DO END DO END DO END DO ! sum of the alerts on all processors IF( lk_mpp ) THEN DO ialert_id = 1, inb_altests CALL mpp_sum(inb_alp(ialert_id)) END DO ENDIF ! print alerts IF( lwp ) THEN ialert_id = 1 ! reference number of this alert cl_alname(ialert_id) = ' NO alerte 1 ' ! name of the alert WRITE(numout,*) ' time step ',kt WRITE(numout,*) ' All alerts at the end of ice model ' DO ialert_id = 1, inb_altests WRITE(numout,*) ialert_id, cl_alname(ialert_id)//' : ', inb_alp(ialert_id), ' times ! ' END DO ENDIF ! END SUBROUTINE lim_ctl SUBROUTINE lim_prt_state( kt, ki, kj, kn, cd1 ) !!----------------------------------------------------------------------- !! *** ROUTINE lim_prt_state *** !! !! ** Purpose : Writes global ice state on the (i,j) point !! in ocean.ouput !! 3 possibilities exist !! n = 1/-1 -> simple ice state (plus Mechanical Check if -1) !! n = 2 -> exhaustive state !! n = 3 -> ice/ocean salt fluxes !! !! ** input : point coordinates (i,j) !! n : number of the option !!------------------------------------------------------------------- INTEGER , INTENT(in) :: kt ! ocean time step INTEGER , INTENT(in) :: ki, kj, kn ! ocean gridpoint indices CHARACTER(len=*), INTENT(in) :: cd1 ! !! INTEGER :: jl, ji, jj !!------------------------------------------------------------------- DO ji = mi0(ki), mi1(ki) DO jj = mj0(kj), mj1(kj) WRITE(numout,*) ' time step ',kt,' ',cd1 ! print title !---------------- ! Simple state !---------------- IF ( kn == 1 .OR. kn == -1 ) THEN WRITE(numout,*) ' lim_prt_state - Point : ',ji,jj WRITE(numout,*) ' ~~~~~~~~~~~~~~ ' WRITE(numout,*) ' Simple state ' WRITE(numout,*) ' masks s,u,v : ', tms(ji,jj), tmu(ji,jj), tmv(ji,jj) WRITE(numout,*) ' lat - long : ', gphit(ji,jj), glamt(ji,jj) WRITE(numout,*) ' Time step : ', numit WRITE(numout,*) ' - Ice drift ' WRITE(numout,*) ' ~~~~~~~~~~~ ' WRITE(numout,*) ' u_ice(i-1,j) : ', u_ice(ji-1,jj) WRITE(numout,*) ' u_ice(i ,j) : ', u_ice(ji,jj) WRITE(numout,*) ' v_ice(i ,j-1): ', v_ice(ji,jj-1) WRITE(numout,*) ' v_ice(i ,j) : ', v_ice(ji,jj) WRITE(numout,*) ' strength : ', strength(ji,jj) WRITE(numout,*) WRITE(numout,*) ' - Cell values ' WRITE(numout,*) ' ~~~~~~~~~~~ ' WRITE(numout,*) ' cell area : ', area(ji,jj) WRITE(numout,*) ' at_i : ', at_i(ji,jj) WRITE(numout,*) ' vt_i : ', vt_i(ji,jj) WRITE(numout,*) ' vt_s : ', vt_s(ji,jj) DO jl = 1, jpl WRITE(numout,*) ' - Category (', jl,')' WRITE(numout,*) ' a_i : ', a_i(ji,jj,jl) WRITE(numout,*) ' ht_i : ', ht_i(ji,jj,jl) WRITE(numout,*) ' ht_s : ', ht_s(ji,jj,jl) WRITE(numout,*) ' v_i : ', v_i(ji,jj,jl) WRITE(numout,*) ' v_s : ', v_s(ji,jj,jl) WRITE(numout,*) ' e_s : ', e_s(ji,jj,1,jl)/1.0e9 WRITE(numout,*) ' e_i : ', e_i(ji,jj,1:nlay_i,jl)/1.0e9 WRITE(numout,*) ' t_su : ', t_su(ji,jj,jl) WRITE(numout,*) ' t_snow : ', t_s(ji,jj,1,jl) WRITE(numout,*) ' t_i : ', t_i(ji,jj,1:nlay_i,jl) WRITE(numout,*) ' sm_i : ', sm_i(ji,jj,jl) WRITE(numout,*) ' smv_i : ', smv_i(ji,jj,jl) WRITE(numout,*) END DO ENDIF IF( kn == -1 ) THEN WRITE(numout,*) ' Mechanical Check ************** ' WRITE(numout,*) ' Check what means ice divergence ' WRITE(numout,*) ' Total ice concentration ', at_i (ji,jj) WRITE(numout,*) ' Total lead fraction ', ato_i(ji,jj) WRITE(numout,*) ' Sum of both ', ato_i(ji,jj) + at_i(ji,jj) WRITE(numout,*) ' Sum of both minus 1 ', ato_i(ji,jj) + at_i(ji,jj) - 1.00 ENDIF !-------------------- ! Exhaustive state !-------------------- IF ( kn .EQ. 2 ) THEN WRITE(numout,*) ' lim_prt_state - Point : ',ji,jj WRITE(numout,*) ' ~~~~~~~~~~~~~~ ' WRITE(numout,*) ' Exhaustive state ' WRITE(numout,*) ' lat - long ', gphit(ji,jj), glamt(ji,jj) WRITE(numout,*) ' Time step ', numit WRITE(numout,*) WRITE(numout,*) ' - Cell values ' WRITE(numout,*) ' ~~~~~~~~~~~ ' WRITE(numout,*) ' cell area : ', area(ji,jj) WRITE(numout,*) ' at_i : ', at_i(ji,jj) WRITE(numout,*) ' vt_i : ', vt_i(ji,jj) WRITE(numout,*) ' vt_s : ', vt_s(ji,jj) WRITE(numout,*) ' u_ice(i-1,j) : ', u_ice(ji-1,jj) WRITE(numout,*) ' u_ice(i ,j) : ', u_ice(ji,jj) WRITE(numout,*) ' v_ice(i ,j-1): ', v_ice(ji,jj-1) WRITE(numout,*) ' v_ice(i ,j) : ', v_ice(ji,jj) WRITE(numout,*) ' strength : ', strength(ji,jj) WRITE(numout,*) ' d_u_ice_dyn : ', d_u_ice_dyn(ji,jj), ' d_v_ice_dyn : ', d_v_ice_dyn(ji,jj) WRITE(numout,*) ' u_ice_b : ', u_ice_b(ji,jj) , ' v_ice_b : ', v_ice_b(ji,jj) WRITE(numout,*) DO jl = 1, jpl WRITE(numout,*) ' - Category (',jl,')' WRITE(numout,*) ' ~~~~~~~~ ' WRITE(numout,*) ' ht_i : ', ht_i(ji,jj,jl) , ' ht_s : ', ht_s(ji,jj,jl) WRITE(numout,*) ' t_i : ', t_i(ji,jj,1:nlay_i,jl) WRITE(numout,*) ' t_su : ', t_su(ji,jj,jl) , ' t_s : ', t_s(ji,jj,1,jl) WRITE(numout,*) ' sm_i : ', sm_i(ji,jj,jl) , ' o_i : ', o_i(ji,jj,jl) WRITE(numout,*) ' a_i : ', a_i(ji,jj,jl) , ' a_i_b : ', a_i_b(ji,jj,jl) WRITE(numout,*) ' d_a_i_trp : ', d_a_i_trp(ji,jj,jl) , ' d_a_i_thd : ', d_a_i_thd(ji,jj,jl) WRITE(numout,*) ' v_i : ', v_i(ji,jj,jl) , ' v_i_b : ', v_i_b(ji,jj,jl) WRITE(numout,*) ' d_v_i_trp : ', d_v_i_trp(ji,jj,jl) , ' d_v_i_thd : ', d_v_i_thd(ji,jj,jl) WRITE(numout,*) ' v_s : ', v_s(ji,jj,jl) , ' v_s_b : ', v_s_b(ji,jj,jl) WRITE(numout,*) ' d_v_s_trp : ', d_v_s_trp(ji,jj,jl) , ' d_v_s_thd : ', d_v_s_thd(ji,jj,jl) WRITE(numout,*) ' e_i1 : ', e_i(ji,jj,1,jl)/1.0e9 , ' ei1 : ', e_i_b(ji,jj,1,jl)/1.0e9 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 WRITE(numout,*) ' e_i2 : ', e_i(ji,jj,2,jl)/1.0e9 , ' ei2_b : ', e_i_b(ji,jj,2,jl)/1.0e9 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 WRITE(numout,*) ' e_snow : ', e_s(ji,jj,1,jl) , ' e_snow_b : ', e_s_b(ji,jj,1,jl) 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) WRITE(numout,*) ' smv_i : ', smv_i(ji,jj,jl) , ' smv_i_b : ', smv_i_b(ji,jj,jl) WRITE(numout,*) ' d_smv_i_trp: ', d_smv_i_trp(ji,jj,jl) , ' d_smv_i_thd: ', d_smv_i_thd(ji,jj,jl) WRITE(numout,*) ' oa_i : ', oa_i(ji,jj,jl) , ' oa_i_b : ', oa_i_b(ji,jj,jl) WRITE(numout,*) ' d_oa_i_trp : ', d_oa_i_trp(ji,jj,jl) , ' d_oa_i_thd : ', d_oa_i_thd(ji,jj,jl) END DO !jl WRITE(numout,*) WRITE(numout,*) ' - Heat / FW fluxes ' WRITE(numout,*) ' ~~~~~~~~~~~~~~~~ ' WRITE(numout,*) ' - Heat fluxes in and out the ice ***' WRITE(numout,*) ' qsr_ini : ', pfrld(ji,jj) * qsr(ji,jj) + SUM( a_i_b(ji,jj,:) * qsr_ice(ji,jj,:) ) WRITE(numout,*) ' qns_ini : ', pfrld(ji,jj) * qns(ji,jj) + SUM( a_i_b(ji,jj,:) * qns_ice(ji,jj,:) ) WRITE(numout,*) WRITE(numout,*) WRITE(numout,*) ' sst : ', sst_m(ji,jj) WRITE(numout,*) ' sss : ', sss_m(ji,jj) WRITE(numout,*) WRITE(numout,*) ' - Stresses ' WRITE(numout,*) ' ~~~~~~~~ ' WRITE(numout,*) ' utau_ice : ', utau_ice(ji,jj) WRITE(numout,*) ' vtau_ice : ', vtau_ice(ji,jj) WRITE(numout,*) ' utau : ', utau (ji,jj) WRITE(numout,*) ' vtau : ', vtau (ji,jj) WRITE(numout,*) ' oc. vel. u : ', u_oce (ji,jj) WRITE(numout,*) ' oc. vel. v : ', v_oce (ji,jj) ENDIF !--------------------- ! Salt / heat fluxes !--------------------- IF ( kn .EQ. 3 ) THEN WRITE(numout,*) ' lim_prt_state - Point : ',ji,jj WRITE(numout,*) ' ~~~~~~~~~~~~~~ ' WRITE(numout,*) ' - Salt / Heat Fluxes ' WRITE(numout,*) ' ~~~~~~~~~~~~~~~~ ' WRITE(numout,*) ' lat - long ', gphit(ji,jj), glamt(ji,jj) WRITE(numout,*) ' Time step ', numit WRITE(numout,*) WRITE(numout,*) ' - Heat fluxes at bottom interface ***' WRITE(numout,*) ' qsr : ', qsr(ji,jj) WRITE(numout,*) ' qns : ', qns(ji,jj) WRITE(numout,*) WRITE(numout,*) ' hfx_mass : ', hfx_thd(ji,jj) + hfx_dyn(ji,jj) + hfx_snw(ji,jj) + hfx_res(ji,jj) WRITE(numout,*) ' hfx_in : ', hfx_in(ji,jj) WRITE(numout,*) ' hfx_out : ', hfx_out(ji,jj) WRITE(numout,*) ' dhc : ', diag_heat_dhc(ji,jj) WRITE(numout,*) WRITE(numout,*) ' hfx_dyn : ', hfx_dyn(ji,jj) WRITE(numout,*) ' hfx_thd : ', hfx_thd(ji,jj) WRITE(numout,*) ' hfx_res : ', hfx_res(ji,jj) WRITE(numout,*) ' fhtur : ', fhtur(ji,jj) WRITE(numout,*) ' qlead : ', qlead(ji,jj) * r1_rdtice WRITE(numout,*) WRITE(numout,*) ' - Salt fluxes at bottom interface ***' WRITE(numout,*) ' emp : ', emp (ji,jj) WRITE(numout,*) ' sfx : ', sfx (ji,jj) WRITE(numout,*) ' sfx_res : ', sfx_res(ji,jj) WRITE(numout,*) ' sfx_bri : ', sfx_bri(ji,jj) WRITE(numout,*) ' sfx_dyn : ', sfx_dyn(ji,jj) WRITE(numout,*) WRITE(numout,*) ' - Momentum fluxes ' WRITE(numout,*) ' utau : ', utau(ji,jj) WRITE(numout,*) ' vtau : ', vtau(ji,jj) ENDIF WRITE(numout,*) ' ' ! END DO END DO ! END SUBROUTINE lim_prt_state FUNCTION fice_cell_ave ( ptab ) !!-------------------------------------------------------------------------- !! * Compute average over categories, for grid cell (ice covered and free ocean) !!-------------------------------------------------------------------------- REAL (wp), DIMENSION (jpi,jpj) :: fice_cell_ave REAL (wp), DIMENSION (jpi,jpj,jpl), INTENT (in) :: ptab INTEGER :: jl ! Dummy loop index fice_cell_ave (:,:) = 0.0_wp DO jl = 1, jpl fice_cell_ave (:,:) = fice_cell_ave (:,:) & & + a_i (:,:,jl) * ptab (:,:,jl) END DO END FUNCTION fice_cell_ave FUNCTION fice_ice_ave ( ptab ) !!-------------------------------------------------------------------------- !! * Compute average over categories, for ice covered part of grid cell !!-------------------------------------------------------------------------- REAL (kind=wp), DIMENSION (jpi,jpj) :: fice_ice_ave REAL (kind=wp), DIMENSION (jpi,jpj,jpl), INTENT(in) :: ptab fice_ice_ave (:,:) = 0.0_wp WHERE ( at_i (:,:) .GT. 0.0_wp ) fice_ice_ave (:,:) = fice_cell_ave ( ptab (:,:,:)) / at_i (:,:) END FUNCTION fice_ice_ave #else !!---------------------------------------------------------------------- !! Default option Dummy module NO LIM 3.0 sea-ice model !!---------------------------------------------------------------------- CONTAINS SUBROUTINE sbc_ice_lim ( kt, kblk ) ! Dummy routine WRITE(*,*) 'sbc_ice_lim: You should not have seen this print! error?', kt, kblk END SUBROUTINE sbc_ice_lim #endif !!====================================================================== END MODULE sbcice_lim