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 !!---------------------------------------------------------------------- #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 !! * 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 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 !!====================================================================== 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) !! INTEGER :: jl ! dummy loop index REAL(wp) :: zcoef ! local scalar REAL(wp), POINTER, DIMENSION(:,:,:) :: zalb_ice_os, zalb_ice_cs ! albedo of the ice under overcast/clear sky REAL(wp), POINTER, DIMENSION(:,:,:) :: zalb_ice ! mean albedo of ice (for coupled) REAL(wp), POINTER, DIMENSION(:,:) :: zalb_ice_all ! Mean albedo over all categories REAL(wp), POINTER, DIMENSION(:,:) :: ztem_ice_all ! Mean temperature over all categories REAL(wp), POINTER, DIMENSION(:,:) :: z_qsr_ice_all ! Mean solar heat flux over all categories REAL(wp), POINTER, DIMENSION(:,:) :: z_qns_ice_all ! Mean non solar heat flux over all categories REAL(wp), POINTER, DIMENSION(:,:) :: z_qla_ice_all ! Mean latent heat flux over all categories REAL(wp), POINTER, DIMENSION(:,:) :: z_dqns_ice_all ! Mean d(qns)/dT over all categories REAL(wp), POINTER, DIMENSION(:,:) :: z_dqla_ice_all ! Mean d(qla)/dT over all categories !!---------------------------------------------------------------------- !- O.M. : why do we allocate all these arrays even when MOD( kt-1, nn_fsbc ) /= 0 ????? IF( nn_timing == 1 ) CALL timing_start('sbc_ice_lim') CALL wrk_alloc( jpi,jpj,jpl, zalb_ice_os, zalb_ice_cs ) #if defined key_coupled IF ( ln_cpl .OR. ln_iceflx_ave .OR. ln_iceflx_linear ) CALL wrk_alloc( jpi,jpj,jpl, zalb_ice) IF ( ln_iceflx_ave .OR. ln_iceflx_linear ) & & CALL wrk_alloc( jpi,jpj, ztem_ice_all, zalb_ice_all, z_qsr_ice_all, z_qns_ice_all, z_qla_ice_all, z_dqns_ice_all, z_dqla_ice_all) #endif 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 = 177 ; jjindx = 112 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(:,:) ! mean surface ocean current at ice velocity point v_oce(:,:) = ssv_m(:,:) ! (C-grid dynamics : U- & V-points as the ocean) ! t_bo(:,:) = tfreez( sss_m ) + rt0 ! masked sea surface freezing temperature [Kelvin] ! ! (set to rt0 over land) CALL albedo_ice( t_su, ht_i, ht_s, zalb_ice_cs, zalb_ice_os ) ! ... ice albedo DO jl = 1, jpl t_su(:,:,jl) = t_su(:,:,jl) + rt0 * ( 1. - tmask(:,:,1) ) END DO IF ( ln_cpl ) zalb_ice (:,:,:) = 0.5 * ( zalb_ice_cs (:,:,:) + zalb_ice_os (:,:,:) ) #if defined key_coupled IF ( ln_iceflx_ave .OR. ln_iceflx_linear ) THEN ! ! Compute mean albedo and temperature zalb_ice_all (:,:) = fice_ice_ave ( zalb_ice (:,:,:) ) ztem_ice_all (:,:) = fice_ice_ave ( tn_ice (:,:,:) ) ! ENDIF #endif ! Bulk formulea - 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( 3 ) ! CLIO bulk formulation CALL blk_ice_clio( t_su , zalb_ice_cs, zalb_ice_os, & & utau_ice , vtau_ice , qns_ice , qsr_ice , & & qla_ice , dqns_ice , dqla_ice , & & tprecip , sprecip , & & fr1_i0 , fr2_i0 , cp_ice_msh, jpl ) ! CASE( 4 ) ! CORE bulk formulation ! MV 2014 ! We must account for cloud fraction in the computation of the albedo ! The present ref just uses the clear sky value ! The overcast sky value is 0.06 higher, and polar skies are mostly overcast ! CORE has no cloud fraction, hence we must prescribe it ! Mean summer cloud fraction computed from CLIO = 0.81 zalb_ice(:,:,:) = 0.19 * zalb_ice_cs(:,:,:) + 0.81 * zalb_ice_os(:,:,:) ! Following line, we replace zalb_ice_cs by simply zalb_ice 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 ) ! CASE ( 5 ) zalb_ice (:,:,:) = 0.5 * ( zalb_ice_cs (:,:,:) + zalb_ice_os (:,:,:) ) CALL sbc_cpl_ice_tau( utau_ice , vtau_ice ) CALL sbc_cpl_ice_flx( p_frld=ato_i, palbi=zalb_ice, psst=sst_m, pist=tn_ice ) ! Latent heat flux is forced to 0 in coupled : ! it is included in qns (non-solar heat flux) qla_ice (:,:,:) = 0.0e0_wp dqla_ice (:,:,:) = 0.0e0_wp ! END SELECT ! Average over all categories #if defined key_coupled IF ( ln_iceflx_ave .OR. ln_iceflx_linear ) THEN z_qns_ice_all (:,:) = fice_ice_ave ( qns_ice (:,:,:) ) z_qsr_ice_all (:,:) = fice_ice_ave ( qsr_ice (:,:,:) ) z_dqns_ice_all (:,:) = fice_ice_ave ( dqns_ice (:,:,:) ) z_qla_ice_all (:,:) = fice_ice_ave ( qla_ice (:,:,:) ) z_dqla_ice_all (:,:) = fice_ice_ave ( dqla_ice (:,:,:) ) DO jl = 1, jpl dqns_ice (:,:,jl) = z_dqns_ice_all (:,:) dqla_ice (:,:,jl) = z_dqla_ice_all (:,:) END DO ! IF ( ln_iceflx_ave ) THEN DO jl = 1, jpl qns_ice (:,:,jl) = z_qns_ice_all (:,:) qsr_ice (:,:,jl) = z_qsr_ice_all (:,:) qla_ice (:,:,jl) = z_qla_ice_all (:,:) END DO END IF ! IF ( ln_iceflx_linear ) THEN DO jl = 1, jpl qns_ice (:,:,jl) = z_qns_ice_all(:,:) + z_dqns_ice_all(:,:) * (tn_ice(:,:,jl) - ztem_ice_all(:,:)) qla_ice (:,:,jl) = z_qla_ice_all(:,:) + z_dqla_ice_all(:,:) * (tn_ice(:,:,jl) - ztem_ice_all(:,:)) qsr_ice (:,:,jl) = (1.0e0_wp-zalb_ice(:,:,jl)) / (1.0e0_wp-zalb_ice_all(:,:)) * z_qsr_ice_all(:,:) END DO END IF END IF #endif ! !----------------------! ! ! LIM-3 time-stepping ! ! !----------------------! ! numit = numit + nn_fsbc ! Ice model time step ! ! ! Store previous ice values !!gm : remark old_... should becomes ...b as tn versus tb old_a_i (:,:,:) = a_i (:,:,:) ! ice area old_e_i (:,:,:,:) = e_i (:,:,:,:) ! ice thermal energy old_v_i (:,:,:) = v_i (:,:,:) ! ice volume old_v_s (:,:,:) = v_s (:,:,:) ! snow volume old_e_s (:,:,:,:) = e_s (:,:,:,:) ! snow thermal energy old_smv_i(:,:,:) = smv_i(:,:,:) ! salt content old_oa_i (:,:,:) = oa_i (:,:,:) ! areal age content ! old_u_ice(:,:) = u_ice(:,:) old_v_ice(:,:) = v_ice(:,:) ! ! intialisation to zero !!gm is it truly necessary ??? d_a_i_thd (:,:,:) = 0._wp ; d_a_i_trp (:,:,:) = 0._wp d_v_i_thd (:,:,:) = 0._wp ; d_v_i_trp (:,:,:) = 0._wp d_e_i_thd (:,:,:,:) = 0._wp ; d_e_i_trp (:,:,:,:) = 0._wp d_v_s_thd (:,:,:) = 0._wp ; d_v_s_trp (:,:,:) = 0._wp d_e_s_thd (:,:,:,:) = 0._wp ; d_e_s_trp (:,:,:,:) = 0._wp d_smv_i_thd(:,:,:) = 0._wp ; d_smv_i_trp(:,:,:) = 0._wp d_oa_i_thd (:,:,:) = 0._wp ; d_oa_i_trp (:,:,:) = 0._wp ! d_u_ice_dyn(:,:) = 0._wp d_v_ice_dyn(:,:) = 0._wp ! sfx (:,:) = 0._wp ; sfx_thd (:,:) = 0._wp sfx_bri(:,:) = 0._wp ; sfx_mec (:,:) = 0._wp ; sfx_res (:,:) = 0._wp fhbri (:,:) = 0._wp ; fheat_mec(:,:) = 0._wp ; fheat_res(:,:) = 0._wp fhmec (:,:) = 0._wp ; fmmec (:,:) = 0._wp fmmflx (:,:) = 0._wp focea2D(:,:) = 0._wp fsup2D (:,:) = 0._wp ! used in limthd.F90 rdvosif(:,:) = 0._wp ! variation of ice volume at surface rdvobif(:,:) = 0._wp ! variation of ice volume at bottom fdvolif(:,:) = 0._wp ! total variation of ice volume rdvonif(:,:) = 0._wp ! lateral variation of ice volume fstric (:,:) = 0._wp ! part of solar radiation transmitted through the ice ffltbif(:,:) = 0._wp ! linked with fstric qfvbq (:,:) = 0._wp ! linked with fstric rdm_snw(:,:) = 0._wp ! variation of snow mass per unit area rdm_ice(:,:) = 0._wp ! variation of ice mass per unit area hicifp (:,:) = 0._wp ! daily thermodynamic ice production. ! diag_sni_gr(:,:) = 0._wp ; diag_lat_gr(:,:) = 0._wp diag_bot_gr(:,:) = 0._wp ; diag_dyn_gr(:,:) = 0._wp diag_bot_me(:,:) = 0._wp ; diag_sur_me(:,:) = 0._wp diag_res_pr(:,:) = 0._wp ; diag_trp_vi(:,:) = 0._wp ! dynamical invariants delta_i(:,:) = 0._wp ; divu_i(:,:) = 0._wp ; shear_i(:,:) = 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 old_u_ice(:,:) = u_ice (:,:) old_v_ice(:,:) = v_ice (:,:) old_a_i(:,:,:) = a_i (:,:,:) old_v_s(:,:,:) = v_s (:,:,:) old_v_i(:,:,:) = v_i (:,:,:) old_e_s(:,:,:,:) = e_s (:,:,:,:) old_e_i(:,:,:,:) = e_i (:,:,:,:) old_oa_i(:,:,:) = oa_i(:,:,:) old_smv_i(:,:,:) = 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(:,:) ! 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 CALL lim_var_glo2eqv ! this CALL is maybe not necessary (Martin) 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 !clem # if ! defined key_iomput CALL lim_wri( 1 ) ! Ice outputs !clem # endif 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 ! 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!!! ! CALL wrk_dealloc( jpi,jpj,jpl, zalb_ice_os, zalb_ice_cs ) #if defined key_coupled IF ( ln_cpl .OR. ln_iceflx_ave .OR. ln_iceflx_linear ) CALL wrk_dealloc( jpi,jpj,jpl, zalb_ice) IF ( ln_iceflx_ave .OR. ln_iceflx_linear ) & & CALL wrk_dealloc( jpi,jpj, ztem_ice_all, zalb_ice_all, z_qsr_ice_all, z_qns_ice_all, z_qla_ice_all, z_dqns_ice_all, z_dqla_ice_all) #endif ! IF( nn_timing == 1 ) CALL timing_stop('sbc_ice_lim') ! END SUBROUTINE sbc_ice_lim 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_old ', a_i (ji,jj,jl), old_a_i (ji,jj,jl) !WRITE(numout,*) ' v_i *** v_i_old ', v_i (ji,jj,jl), old_v_i (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,*) ' s_i_newice : ', s_i_newice(ji,jj,1:jpl) ! 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) , ' old_a_i : ', old_a_i (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) , ' old_v_i : ', old_v_i (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) !WRITE(numout,*) ' qcmif : ', qcmif(ji,jj) !WRITE(numout,*) ' qldif : ', qldif(ji,jj) !WRITE(numout,*) ' qcmif : ', qcmif(ji,jj) / rdt_ice !WRITE(numout,*) ' qldif : ', qldif(ji,jj) / rdt_ice !WRITE(numout,*) ' qfvbq : ', qfvbq(ji,jj) !WRITE(numout,*) ' qdtcn : ', qdtcn(ji,jj) !WRITE(numout,*) ' qfvbq / dt: ', qfvbq(ji,jj) / rdt_ice !WRITE(numout,*) ' qdtcn / dt: ', qdtcn(ji,jj) / rdt_ice !WRITE(numout,*) ' fdtcn : ', fdtcn(ji,jj) !WRITE(numout,*) ' fhmec : ', fhmec(ji,jj) !WRITE(numout,*) ' fheat_mec : ', fheat_mec(ji,jj) !WRITE(numout,*) ' fheat_res : ', fheat_res(ji,jj) !WRITE(numout,*) ' fhbri : ', fhbri(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,*) ' old_u_ice : ', old_u_ice(ji,jj) , ' old_v_ice : ', old_v_ice(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) , ' old_a_i : ', old_a_i(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) , ' old_v_i : ', old_v_i(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) , ' old_v_s : ', old_v_s(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 , ' old_ei1 : ', old_e_i(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 , ' old_ei2 : ', old_e_i(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) , ' old_e_snow : ', old_e_s(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) , ' old_smv_i : ', old_smv_i(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) , ' old_oa_i : ', old_oa_i(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,*) ' emp : ', emp (ji,jj) WRITE(numout,*) ' sfx : ', sfx (ji,jj) WRITE(numout,*) ' sfx_thd : ', sfx_thd(ji,jj) WRITE(numout,*) ' sfx_bri : ', sfx_bri (ji,jj) WRITE(numout,*) ' sfx_mec : ', sfx_mec (ji,jj) WRITE(numout,*) ' sfx_res : ', sfx_res(ji,jj) WRITE(numout,*) ' fmmec : ', fmmec (ji,jj) WRITE(numout,*) ' fhmec : ', fhmec (ji,jj) WRITE(numout,*) ' fhbri : ', fhbri (ji,jj) WRITE(numout,*) ' fheat_mec : ', fheat_mec(ji,jj) 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,*) ' fdtcn : ', fdtcn(ji,jj) WRITE(numout,*) ' qcmif : ', qcmif(ji,jj) * r1_rdtice WRITE(numout,*) ' qldif : ', qldif(ji,jj) * r1_rdtice WRITE(numout,*) WRITE(numout,*) ' - Salt fluxes at bottom interface ***' WRITE(numout,*) ' emp : ', emp (ji,jj) WRITE(numout,*) ' sfx_bri : ', sfx_bri(ji,jj) WRITE(numout,*) ' sfx : ', sfx (ji,jj) WRITE(numout,*) ' sfx_res : ', sfx_res(ji,jj) WRITE(numout,*) ' sfx_mec : ', sfx_mec(ji,jj) WRITE(numout,*) ' - Heat fluxes at bottom interface ***' WRITE(numout,*) ' fheat_res : ', fheat_res(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 #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