MODULE limtrp !!====================================================================== !! *** MODULE limtrp *** !! LIM transport ice model : sea-ice advection/diffusion !!====================================================================== !! History : LIM-2 ! 2000-01 (M.A. Morales Maqueda, H. Goosse, and T. Fichefet) Original code !! 3.0 ! 2005-11 (M. Vancoppenolle) Multi-layer sea ice, salinity variations !! 4.0 ! 2011-02 (G. Madec) dynamical allocation !!---------------------------------------------------------------------- #if defined key_lim3 !!---------------------------------------------------------------------- !! 'key_lim3' LIM3 sea-ice model !!---------------------------------------------------------------------- !! lim_trp : advection/diffusion process of sea ice !!---------------------------------------------------------------------- USE phycst ! physical constant USE dom_oce ! ocean domain USE sbc_oce ! ocean surface boundary condition USE par_ice ! ice parameter USE dom_ice ! ice domain USE ice ! ice variables USE limadv ! ice advection USE limhdf ! ice horizontal diffusion USE in_out_manager ! I/O manager USE lbclnk ! lateral boundary conditions -- MPP exchanges USE lib_mpp ! MPP library USE wrk_nemo ! work arrays USE prtctl ! Print control USE lib_fortran ! Fortran utilities (allows no signed zero when 'key_nosignedzero' defined) USE limvar ! clem for ice thickness correction USE timing ! Timing USE limcons ! conservation tests IMPLICIT NONE PRIVATE PUBLIC lim_trp ! called by ice_step !! * Substitution # include "vectopt_loop_substitute.h90" !!---------------------------------------------------------------------- !! NEMO/LIM3 4.0 , UCL - NEMO Consortium (2011) !! $Id$ !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) !!---------------------------------------------------------------------- CONTAINS SUBROUTINE lim_trp( kt ) !!------------------------------------------------------------------- !! *** ROUTINE lim_trp *** !! !! ** purpose : advection/diffusion process of sea ice !! !! ** method : variables included in the process are scalar, !! other values are considered as second order. !! For advection, a second order Prather scheme is used. !! !! ** action : !!--------------------------------------------------------------------- INTEGER, INTENT(in) :: kt ! number of iteration ! INTEGER :: ji, jj, jk, jl, jn ! dummy loop indices INTEGER :: initad ! number of sub-timestep for the advection REAL(wp) :: zcfl , zusnit ! - - ! REAL(wp), POINTER, DIMENSION(:,:) :: zui_u, zvi_v, zsm, zs0at, zs0ow REAL(wp), POINTER, DIMENSION(:,:,:) :: zs0ice, zs0sn, zs0a, zs0c0 , zs0sm , zs0oi REAL(wp), POINTER, DIMENSION(:,:,:,:) :: zs0e REAL(wp), POINTER, DIMENSION(:,:,:) :: zviold, zvsold ! old ice volume... REAL(wp), POINTER, DIMENSION(:,:,:) :: zaiold, zhimax ! old ice concentration and thickness REAL(wp), POINTER, DIMENSION(:,:) :: zeiold, zesold ! old enthalpies REAL(wp) :: zdv, zda, zvi, zvs, zsmv, zes, zei ! REAL(wp) :: zvi_b, zsmv_b, zei_b, zfs_b, zfw_b, zft_b !!--------------------------------------------------------------------- IF( nn_timing == 1 ) CALL timing_start('limtrp') CALL wrk_alloc( jpi, jpj, zui_u, zvi_v, zsm, zs0at, zs0ow, zeiold, zesold ) CALL wrk_alloc( jpi, jpj, jpl, zs0ice, zs0sn, zs0a, zs0c0 , zs0sm , zs0oi ) CALL wrk_alloc( jpi, jpj, nlay_i+1, jpl, zs0e ) CALL wrk_alloc( jpi, jpj, jpl, zaiold, zhimax, zviold, zvsold ) ! clem IF( numit == nstart .AND. lwp ) THEN WRITE(numout,*) IF( ln_limdyn ) THEN ; WRITE(numout,*) 'lim_trp : Ice transport ' ELSE ; WRITE(numout,*) 'lim_trp : No ice advection as ln_limdyn = ', ln_limdyn ENDIF WRITE(numout,*) '~~~~~~~~~~~~' ENDIF zsm(:,:) = area(:,:) ! !-------------------------------------! IF( ln_limdyn ) THEN ! Advection of sea ice properties ! ! !-------------------------------------! ! conservation test IF( ln_limdiahsb ) CALL lim_cons_hsm(0, 'limtrp', zvi_b, zsmv_b, zei_b, zfw_b, zfs_b, zft_b) ! mass and salt flux init (clem) zviold(:,:,:) = v_i(:,:,:) zeiold(:,:) = SUM( SUM( e_i(:,:,1:nlay_i,:), dim=4 ), dim=3 ) zesold(:,:) = SUM( SUM( e_s(:,:,1:nlay_s,:), dim=4 ), dim=3 ) !--- Thickness correction init. (clem) ------------------------------- CALL lim_var_glo2eqv zaiold(:,:,:) = a_i(:,:,:) !--------------------------------------------------------------------- ! Record max of the surrounding ice thicknesses for correction in limupdate ! in case advection creates ice too thick. !--------------------------------------------------------------------- zhimax(:,:,:) = ht_i(:,:,:) DO jl = 1, jpl DO jj = 2, jpjm1 DO ji = 2, jpim1 zhimax(ji,jj,jl) = MAXVAL( ht_i(ji-1:ji+1,jj-1:jj+1,jl) ) !zhimax(ji,jj,jl) = ( ht_i(ji ,jj ,jl) * tmask(ji, jj ,1) + ht_i(ji-1,jj-1,jl) * tmask(ji-1,jj-1,1) + ht_i(ji+1,jj+1,jl) * tmask(ji+1,jj+1,1) & ! & + ht_i(ji-1,jj ,jl) * tmask(ji-1,jj ,1) + ht_i(ji ,jj-1,jl) * tmask(ji ,jj-1,1) & ! & + ht_i(ji+1,jj ,jl) * tmask(ji+1,jj ,1) + ht_i(ji ,jj+1,jl) * tmask(ji ,jj+1,1) & ! & + ht_i(ji-1,jj+1,jl) * tmask(ji-1,jj+1,1) + ht_i(ji+1,jj-1,jl) * tmask(ji+1,jj-1,1) ) END DO END DO CALL lbc_lnk(zhimax(:,:,jl),'T',1.) END DO !------------------------- ! transported fields !------------------------- ! Snow vol, ice vol, salt and age contents, area zs0ow(:,:) = ato_i(:,:) * area(:,:) ! Open water area DO jl = 1, jpl zs0sn (:,:,jl) = v_s (:,:,jl) * area(:,:) ! Snow volume zs0ice(:,:,jl) = v_i (:,:,jl) * area(:,:) ! Ice volume zs0a (:,:,jl) = a_i (:,:,jl) * area(:,:) ! Ice area zs0sm (:,:,jl) = smv_i(:,:,jl) * area(:,:) ! Salt content zs0oi (:,:,jl) = oa_i (:,:,jl) * area(:,:) ! Age content zs0c0 (:,:,jl) = e_s (:,:,1,jl) ! Snow heat content zs0e (:,:,:,jl) = e_i (:,:,:,jl) ! Ice heat content END DO !-------------------------- ! Advection of Ice fields (Prather scheme) !-------------------------- ! If ice drift field is too fast, use an appropriate time step for advection. ! CFL test for stability zcfl = MAXVAL( ABS( u_ice(:,:) ) * rdt_ice / e1u(:,:) ) zcfl = MAX( zcfl, MAXVAL( ABS( v_ice(:,:) ) * rdt_ice / e2v(:,:) ) ) IF(lk_mpp ) CALL mpp_max( zcfl ) !!gm more readability: ! IF( zcfl > 0.5 ) THEN ; initad = 2 ; zusnit = 0.5_wp ! ELSE ; initad = 1 ; zusnit = 1.0_wp ! ENDIF !!gm end initad = 1 + NINT( MAX( 0._wp, SIGN( 1._wp, zcfl-0.5 ) ) ) zusnit = 1.0 / REAL( initad ) IF( zcfl > 0.5 .AND. lwp ) & WRITE(numout,*) 'lim_trp : CFL violation at day ', nday, ', cfl = ', zcfl, & & ': the ice time stepping is split in two' IF( MOD( ( kt - 1) / nn_fsbc , 2 ) == 0 ) THEN !== odd ice time step: adv_x then adv_y ==! DO jn = 1,initad CALL lim_adv_x( zusnit, u_ice, 1._wp , zsm, zs0ow (:,:), sxopw(:,:), & !--- ice open water area & sxxopw(:,:), syopw(:,:), syyopw(:,:), sxyopw(:,:) ) CALL lim_adv_y( zusnit, v_ice, 0._wp, zsm, zs0ow (:,:), sxopw(:,:), & & sxxopw(:,:), syopw(:,:), syyopw(:,:), sxyopw(:,:) ) DO jl = 1, jpl CALL lim_adv_x( zusnit, u_ice, 1._wp , zsm, zs0ice(:,:,jl), sxice(:,:,jl), & !--- ice volume --- & sxxice(:,:,jl), syice(:,:,jl), syyice(:,:,jl), sxyice(:,:,jl) ) CALL lim_adv_y( zusnit, v_ice, 0._wp, zsm, zs0ice(:,:,jl), sxice(:,:,jl), & & sxxice(:,:,jl), syice(:,:,jl), syyice(:,:,jl), sxyice(:,:,jl) ) CALL lim_adv_x( zusnit, u_ice, 1._wp , zsm, zs0sn (:,:,jl), sxsn (:,:,jl), & !--- snow volume --- & sxxsn (:,:,jl), sysn (:,:,jl), syysn (:,:,jl), sxysn (:,:,jl) ) CALL lim_adv_y( zusnit, v_ice, 0._wp, zsm, zs0sn (:,:,jl), sxsn (:,:,jl), & & sxxsn (:,:,jl), sysn (:,:,jl), syysn (:,:,jl), sxysn (:,:,jl) ) CALL lim_adv_x( zusnit, u_ice, 1._wp , zsm, zs0sm (:,:,jl), sxsal(:,:,jl), & !--- ice salinity --- & sxxsal(:,:,jl), sysal(:,:,jl), syysal(:,:,jl), sxysal(:,:,jl) ) CALL lim_adv_y( zusnit, v_ice, 0._wp, zsm, zs0sm (:,:,jl), sxsal(:,:,jl), & & sxxsal(:,:,jl), sysal(:,:,jl), syysal(:,:,jl), sxysal(:,:,jl) ) CALL lim_adv_x( zusnit, u_ice, 1._wp , zsm, zs0oi (:,:,jl), sxage(:,:,jl), & !--- ice age --- & sxxage(:,:,jl), syage(:,:,jl), syyage(:,:,jl), sxyage(:,:,jl) ) CALL lim_adv_y( zusnit, v_ice, 0._wp, zsm, zs0oi (:,:,jl), sxage(:,:,jl), & & sxxage(:,:,jl), syage(:,:,jl), syyage(:,:,jl), sxyage(:,:,jl) ) CALL lim_adv_x( zusnit, u_ice, 1._wp , zsm, zs0a (:,:,jl), sxa (:,:,jl), & !--- ice concentrations --- & sxxa (:,:,jl), sya (:,:,jl), syya (:,:,jl), sxya (:,:,jl) ) CALL lim_adv_y( zusnit, v_ice, 0._wp, zsm, zs0a (:,:,jl), sxa (:,:,jl), & & sxxa (:,:,jl), sya (:,:,jl), syya (:,:,jl), sxya (:,:,jl) ) CALL lim_adv_x( zusnit, u_ice, 1._wp , zsm, zs0c0 (:,:,jl), sxc0 (:,:,jl), & !--- snow heat contents --- & sxxc0 (:,:,jl), syc0 (:,:,jl), syyc0 (:,:,jl), sxyc0 (:,:,jl) ) CALL lim_adv_y( zusnit, v_ice, 0._wp, zsm, zs0c0 (:,:,jl), sxc0 (:,:,jl), & & sxxc0 (:,:,jl), syc0 (:,:,jl), syyc0 (:,:,jl), sxyc0 (:,:,jl) ) DO jk = 1, nlay_i !--- ice heat contents --- CALL lim_adv_x( zusnit, u_ice, 1._wp , zsm, zs0e(:,:,jk,jl), sxe (:,:,jk,jl), & & sxxe(:,:,jk,jl), sye (:,:,jk,jl), & & syye(:,:,jk,jl), sxye(:,:,jk,jl) ) CALL lim_adv_y( zusnit, v_ice, 0._wp, zsm, zs0e(:,:,jk,jl), sxe (:,:,jk,jl), & & sxxe(:,:,jk,jl), sye (:,:,jk,jl), & & syye(:,:,jk,jl), sxye(:,:,jk,jl) ) END DO END DO END DO ELSE DO jn = 1, initad CALL lim_adv_y( zusnit, v_ice, 1._wp , zsm, zs0ow (:,:), sxopw(:,:), & !--- ice open water area & sxxopw(:,:), syopw(:,:), syyopw(:,:), sxyopw(:,:) ) CALL lim_adv_x( zusnit, u_ice, 0._wp, zsm, zs0ow (:,:), sxopw(:,:), & & sxxopw(:,:), syopw(:,:), syyopw(:,:), sxyopw(:,:) ) DO jl = 1, jpl CALL lim_adv_y( zusnit, v_ice, 1._wp , zsm, zs0ice(:,:,jl), sxice(:,:,jl), & !--- ice volume --- & sxxice(:,:,jl), syice(:,:,jl), syyice(:,:,jl), sxyice(:,:,jl) ) CALL lim_adv_x( zusnit, u_ice, 0._wp, zsm, zs0ice(:,:,jl), sxice(:,:,jl), & & sxxice(:,:,jl), syice(:,:,jl), syyice(:,:,jl), sxyice(:,:,jl) ) CALL lim_adv_y( zusnit, v_ice, 1._wp , zsm, zs0sn (:,:,jl), sxsn (:,:,jl), & !--- snow volume --- & sxxsn (:,:,jl), sysn (:,:,jl), syysn (:,:,jl), sxysn (:,:,jl) ) CALL lim_adv_x( zusnit, u_ice, 0._wp, zsm, zs0sn (:,:,jl), sxsn (:,:,jl), & & sxxsn (:,:,jl), sysn (:,:,jl), syysn (:,:,jl), sxysn (:,:,jl) ) CALL lim_adv_y( zusnit, v_ice, 1._wp , zsm, zs0sm (:,:,jl), sxsal(:,:,jl), & !--- ice salinity --- & sxxsal(:,:,jl), sysal(:,:,jl), syysal(:,:,jl), sxysal(:,:,jl) ) CALL lim_adv_x( zusnit, u_ice, 0._wp, zsm, zs0sm (:,:,jl), sxsal(:,:,jl), & & sxxsal(:,:,jl), sysal(:,:,jl), syysal(:,:,jl), sxysal(:,:,jl) ) CALL lim_adv_y( zusnit, v_ice, 1._wp , zsm, zs0oi (:,:,jl), sxage(:,:,jl), & !--- ice age --- & sxxage(:,:,jl), syage(:,:,jl), syyage(:,:,jl), sxyage(:,:,jl) ) CALL lim_adv_x( zusnit, u_ice, 0._wp, zsm, zs0oi (:,:,jl), sxage(:,:,jl), & & sxxage(:,:,jl), syage(:,:,jl), syyage(:,:,jl), sxyage(:,:,jl) ) CALL lim_adv_y( zusnit, v_ice, 1._wp , zsm, zs0a (:,:,jl), sxa (:,:,jl), & !--- ice concentrations --- & sxxa (:,:,jl), sya (:,:,jl), syya (:,:,jl), sxya (:,:,jl) ) CALL lim_adv_x( zusnit, u_ice, 0._wp, zsm, zs0a (:,:,jl), sxa (:,:,jl), & & sxxa (:,:,jl), sya (:,:,jl), syya (:,:,jl), sxya (:,:,jl) ) CALL lim_adv_y( zusnit, v_ice, 1._wp , zsm, zs0c0 (:,:,jl), sxc0 (:,:,jl), & !--- snow heat contents --- & sxxc0 (:,:,jl), syc0 (:,:,jl), syyc0 (:,:,jl), sxyc0 (:,:,jl) ) CALL lim_adv_x( zusnit, u_ice, 0._wp, zsm, zs0c0 (:,:,jl), sxc0 (:,:,jl), & & sxxc0 (:,:,jl), syc0 (:,:,jl), syyc0 (:,:,jl), sxyc0 (:,:,jl) ) DO jk = 1, nlay_i !--- ice heat contents --- CALL lim_adv_y( zusnit, v_ice, 1._wp , zsm, zs0e(:,:,jk,jl), sxe (:,:,jk,jl), & & sxxe(:,:,jk,jl), sye (:,:,jk,jl), & & syye(:,:,jk,jl), sxye(:,:,jk,jl) ) CALL lim_adv_x( zusnit, u_ice, 0._wp, zsm, zs0e(:,:,jk,jl), sxe (:,:,jk,jl), & & sxxe(:,:,jk,jl), sye (:,:,jk,jl), & & syye(:,:,jk,jl), sxye(:,:,jk,jl) ) END DO END DO END DO ENDIF !------------------------------------------- ! Recover the properties from their contents !------------------------------------------- zs0ow(:,:) = zs0ow(:,:) / area(:,:) DO jl = 1, jpl zs0ice(:,:,jl) = zs0ice(:,:,jl) / area(:,:) zs0sn (:,:,jl) = zs0sn (:,:,jl) / area(:,:) zs0sm (:,:,jl) = zs0sm (:,:,jl) / area(:,:) zs0oi (:,:,jl) = zs0oi (:,:,jl) / area(:,:) zs0a (:,:,jl) = zs0a (:,:,jl) / area(:,:) ! END DO !------------------------------------------------------------------------------! ! 4) Diffusion of Ice fields !------------------------------------------------------------------------------! !-------------------------------- ! diffusion of open water area !-------------------------------- zs0at(:,:) = zs0a(:,:,1) ! total ice fraction DO jl = 2, jpl zs0at(:,:) = zs0at(:,:) + zs0a(:,:,jl) END DO ! ! ! Masked eddy diffusivity coefficient at ocean U- and V-points DO jj = 1, jpjm1 ! NB: has not to be defined on jpj line and jpi row DO ji = 1 , fs_jpim1 ! vector opt. pahu(ji,jj) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -zs0at(ji ,jj) ) ) ) & & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -zs0at(ji+1,jj) ) ) ) * ahiu(ji,jj) pahv(ji,jj) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -zs0at(ji,jj ) ) ) ) & & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp,- zs0at(ji,jj+1) ) ) ) * ahiv(ji,jj) END DO END DO ! CALL lim_hdf( zs0ow (:,:) ) ! Diffusion !------------------------------------ ! Diffusion of other ice variables !------------------------------------ DO jl = 1, jpl ! ! Masked eddy diffusivity coefficient at ocean U- and V-points DO jj = 1, jpjm1 ! NB: has not to be defined on jpj line and jpi row DO ji = 1 , fs_jpim1 ! vector opt. pahu(ji,jj) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -zs0a(ji ,jj,jl) ) ) ) & & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -zs0a(ji+1,jj,jl) ) ) ) * ahiu(ji,jj) pahv(ji,jj) = ( 1._wp - MAX( 0._wp, SIGN( 1._wp, -zs0a(ji,jj ,jl) ) ) ) & & * ( 1._wp - MAX( 0._wp, SIGN( 1._wp,- zs0a(ji,jj+1,jl) ) ) ) * ahiv(ji,jj) END DO END DO CALL lim_hdf( zs0ice (:,:,jl) ) CALL lim_hdf( zs0sn (:,:,jl) ) CALL lim_hdf( zs0sm (:,:,jl) ) CALL lim_hdf( zs0oi (:,:,jl) ) CALL lim_hdf( zs0a (:,:,jl) ) CALL lim_hdf( zs0c0 (:,:,jl) ) DO jk = 1, nlay_i CALL lim_hdf( zs0e (:,:,jk,jl) ) END DO END DO !------------------------------------------------------------------------------! ! 5) Update and limit ice properties after transport !------------------------------------------------------------------------------! !-------------------------------------------------- ! 5.1) Recover mean values over the grid squares. !-------------------------------------------------- zs0at(:,:) = 0._wp DO jl = 1, jpl DO jj = 1, jpj DO ji = 1, jpi zs0sn (ji,jj,jl) = MAX( 0._wp, zs0sn (ji,jj,jl) ) zs0ice(ji,jj,jl) = MAX( 0._wp, zs0ice(ji,jj,jl) ) zs0sm (ji,jj,jl) = MAX( 0._wp, zs0sm (ji,jj,jl) ) zs0oi (ji,jj,jl) = MAX( 0._wp, zs0oi (ji,jj,jl) ) zs0a (ji,jj,jl) = MAX( 0._wp, zs0a (ji,jj,jl) ) zs0c0 (ji,jj,jl) = MAX( 0._wp, zs0c0 (ji,jj,jl) ) zs0at (ji,jj) = zs0at(ji,jj) + zs0a(ji,jj,jl) END DO END DO END DO !--------------------------------------------------------- ! 5.2) Update and mask variables !--------------------------------------------------------- DO jl = 1, jpl DO jj = 1, jpj DO ji = 1, jpi rswitch = MAX( 0._wp , SIGN( 1._wp, zs0a(ji,jj,jl) - epsi10 ) ) zvi = zs0ice(ji,jj,jl) zvs = zs0sn (ji,jj,jl) zes = zs0c0 (ji,jj,jl) zsmv = zs0sm (ji,jj,jl) ! ! Remove very small areas v_s(ji,jj,jl) = rswitch * zs0sn (ji,jj,jl) v_i(ji,jj,jl) = rswitch * zs0ice(ji,jj,jl) a_i(ji,jj,jl) = rswitch * zs0a (ji,jj,jl) e_s(ji,jj,1,jl) = rswitch * zs0c0 (ji,jj,jl) ! Ice salinity and age IF( num_sal == 2 ) THEN smv_i(ji,jj,jl) = MAX( MIN( s_i_max * v_i(ji,jj,jl), zsmv ), s_i_min * v_i(ji,jj,jl) ) ENDIF oa_i(ji,jj,jl) = MAX( rswitch * zs0oi(ji,jj,jl) / MAX( a_i(ji,jj,jl), epsi10 ), 0._wp ) * a_i(ji,jj,jl) ! Update fluxes wfx_res(ji,jj) = wfx_res(ji,jj) - ( v_i(ji,jj,jl) - zvi ) * rhoic * r1_rdtice wfx_snw(ji,jj) = wfx_snw(ji,jj) - ( v_s(ji,jj,jl) - zvs ) * rhosn * r1_rdtice sfx_res(ji,jj) = sfx_res(ji,jj) - ( smv_i(ji,jj,jl) - zsmv ) * rhoic * r1_rdtice hfx_res(ji,jj) = hfx_res(ji,jj) + ( e_s(ji,jj,1,jl) - zes ) * unit_fac / area(ji,jj) * r1_rdtice ! W.m-2 <0 END DO END DO END DO DO jl = 1, jpl DO jk = 1, nlay_i DO jj = 1, jpj DO ji = 1, jpi rswitch = MAX( 0._wp , SIGN( 1._wp, zs0a(ji,jj,jl) - epsi10 ) ) zei = zs0e(ji,jj,jk,jl) e_i(ji,jj,jk,jl) = rswitch * MAX( 0._wp, zs0e(ji,jj,jk,jl) ) ! Update fluxes hfx_res(ji,jj) = hfx_res(ji,jj) + ( e_i(ji,jj,jk,jl) - zei ) * unit_fac / area(ji,jj) * r1_rdtice ! W.m-2 <0 END DO !ji END DO ! jj END DO ! jk END DO ! jl !--- Thickness correction in case too high (clem) -------------------------------------------------------- CALL lim_var_glo2eqv DO jl = 1, jpl DO jj = 1, jpj DO ji = 1, jpi IF ( v_i(ji,jj,jl) > 0._wp ) THEN zvi = v_i (ji,jj,jl) zvs = v_s (ji,jj,jl) zsmv = smv_i(ji,jj,jl) zes = e_s (ji,jj,1,jl) zei = SUM( e_i(ji,jj,1:nlay_i,jl) ) zdv = v_i(ji,jj,jl) - zviold(ji,jj,jl) !zda = a_i(ji,jj,jl) - zaiold(ji,jj,jl) rswitch = 1._wp IF ( ( zdv > 0.0 .AND. ht_i(ji,jj,jl) > zhimax(ji,jj,jl) .AND. SUM( zaiold(ji,jj,1:jpl) ) < 0.80 ) .OR. & & ( zdv < 0.0 .AND. ht_i(ji,jj,jl) > zhimax(ji,jj,jl) ) ) THEN ht_i(ji,jj,jl) = MIN( zhimax(ji,jj,jl), hi_max(jl) ) rswitch = MAX( 0._wp, SIGN( 1._wp, ht_i(ji,jj,jl) - epsi20 ) ) a_i(ji,jj,jl) = rswitch * v_i(ji,jj,jl) / MAX( ht_i(ji,jj,jl), epsi20 ) ELSE ht_i(ji,jj,jl) = MAX( MIN( ht_i(ji,jj,jl), hi_max(jl) ), hi_max(jl-1) ) rswitch = MAX( 0._wp, SIGN( 1._wp, ht_i(ji,jj,jl) - epsi20 ) ) a_i(ji,jj,jl) = rswitch * v_i(ji,jj,jl) / MAX( ht_i(ji,jj,jl), epsi20 ) ENDIF ! small correction due to *rswitch for a_i v_i (ji,jj,jl) = rswitch * v_i (ji,jj,jl) v_s (ji,jj,jl) = rswitch * v_s (ji,jj,jl) smv_i(ji,jj,jl) = rswitch * smv_i(ji,jj,jl) e_s(ji,jj,1,jl) = rswitch * e_s(ji,jj,1,jl) e_i(ji,jj,1:nlay_i,jl) = rswitch * e_i(ji,jj,1:nlay_i,jl) ! Update mass fluxes wfx_res(ji,jj) = wfx_res(ji,jj) - ( v_i(ji,jj,jl) - zvi ) * rhoic * r1_rdtice wfx_snw(ji,jj) = wfx_snw(ji,jj) - ( v_s(ji,jj,jl) - zvs ) * rhosn * r1_rdtice sfx_res(ji,jj) = sfx_res(ji,jj) - ( smv_i(ji,jj,jl) - zsmv ) * rhoic * r1_rdtice hfx_res(ji,jj) = hfx_res(ji,jj) + ( e_s(ji,jj,1,jl) - zes ) * unit_fac / area(ji,jj) * r1_rdtice ! W.m-2 <0 hfx_res(ji,jj) = hfx_res(ji,jj) + ( SUM( e_i(ji,jj,1:nlay_i,jl) ) - zei ) * unit_fac / area(ji,jj) * r1_rdtice ! W.m-2 <0 ENDIF END DO END DO END DO ! ------------------------------------------------- !-------------------------------------- ! Impose a_i < amax in mono-category !-------------------------------------- ! IF ( ( nn_monocat == 2 ) .AND. ( jpl == 1 ) ) THEN ! simple conservative piling, comparable with LIM2 DO jj = 1, jpj DO ji = 1, jpi a_i(ji,jj,1) = MIN( a_i(ji,jj,1), amax ) END DO END DO ENDIF ! --- diags --- DO jj = 1, jpj DO ji = 1, jpi diag_trp_ei(ji,jj) = ( SUM( e_i(ji,jj,1:nlay_i,:) ) - zeiold(ji,jj) ) / area(ji,jj) * unit_fac * r1_rdtice diag_trp_es(ji,jj) = ( SUM( e_s(ji,jj,1:nlay_s,:) ) - zesold(ji,jj) ) / area(ji,jj) * unit_fac * r1_rdtice diag_trp_vi(ji,jj) = SUM( v_i(ji,jj,:) - zviold(ji,jj,:) ) * r1_rdtice diag_trp_vs(ji,jj) = SUM( v_s(ji,jj,:) - zvsold(ji,jj,:) ) * r1_rdtice END DO END DO ! --- agglomerate variables ----------------- vt_i (:,:) = 0._wp vt_s (:,:) = 0._wp at_i (:,:) = 0._wp ! DO jl = 1, jpl DO jj = 1, jpj DO ji = 1, jpi ! vt_i(ji,jj) = vt_i(ji,jj) + v_i(ji,jj,jl) ! ice volume vt_s(ji,jj) = vt_s(ji,jj) + v_s(ji,jj,jl) ! snow volume at_i(ji,jj) = at_i(ji,jj) + a_i(ji,jj,jl) ! ice concentration END DO END DO END DO ! ------------------------------------------------- ! open water DO jj = 1, jpj DO ji = 1, jpi ! open water = 1 if at_i=0 rswitch = MAX( 0._wp , SIGN( 1._wp, - at_i(ji,jj) ) ) ato_i(ji,jj) = rswitch + (1._wp - rswitch ) * zs0ow(ji,jj) END DO END DO ! conservation test IF( ln_limdiahsb ) CALL lim_cons_hsm(1, 'limtrp', zvi_b, zsmv_b, zei_b, zfw_b, zfs_b, zft_b) ENDIF IF(ln_ctl) THEN ! Control print CALL prt_ctl_info(' ') CALL prt_ctl_info(' - Cell values : ') CALL prt_ctl_info(' ~~~~~~~~~~~~~ ') CALL prt_ctl(tab2d_1=area , clinfo1=' lim_trp : cell area :') CALL prt_ctl(tab2d_1=at_i , clinfo1=' lim_trp : at_i :') CALL prt_ctl(tab2d_1=vt_i , clinfo1=' lim_trp : vt_i :') CALL prt_ctl(tab2d_1=vt_s , clinfo1=' lim_trp : vt_s :') DO jl = 1, jpl CALL prt_ctl_info(' ') CALL prt_ctl_info(' - Category : ', ivar1=jl) CALL prt_ctl_info(' ~~~~~~~~~~') CALL prt_ctl(tab2d_1=a_i (:,:,jl) , clinfo1= ' lim_trp : a_i : ') CALL prt_ctl(tab2d_1=ht_i (:,:,jl) , clinfo1= ' lim_trp : ht_i : ') CALL prt_ctl(tab2d_1=ht_s (:,:,jl) , clinfo1= ' lim_trp : ht_s : ') CALL prt_ctl(tab2d_1=v_i (:,:,jl) , clinfo1= ' lim_trp : v_i : ') CALL prt_ctl(tab2d_1=v_s (:,:,jl) , clinfo1= ' lim_trp : v_s : ') CALL prt_ctl(tab2d_1=e_s (:,:,1,jl) , clinfo1= ' lim_trp : e_s : ') CALL prt_ctl(tab2d_1=t_su (:,:,jl) , clinfo1= ' lim_trp : t_su : ') CALL prt_ctl(tab2d_1=t_s (:,:,1,jl) , clinfo1= ' lim_trp : t_snow : ') CALL prt_ctl(tab2d_1=sm_i (:,:,jl) , clinfo1= ' lim_trp : sm_i : ') CALL prt_ctl(tab2d_1=smv_i (:,:,jl) , clinfo1= ' lim_trp : smv_i : ') DO jk = 1, nlay_i CALL prt_ctl_info(' ') CALL prt_ctl_info(' - Layer : ', ivar1=jk) CALL prt_ctl_info(' ~~~~~~~') CALL prt_ctl(tab2d_1=t_i(:,:,jk,jl) , clinfo1= ' lim_trp : t_i : ') CALL prt_ctl(tab2d_1=e_i(:,:,jk,jl) , clinfo1= ' lim_trp : e_i : ') END DO END DO ENDIF ! CALL wrk_dealloc( jpi, jpj, zui_u, zvi_v, zsm, zs0at, zs0ow, zeiold, zesold ) CALL wrk_dealloc( jpi, jpj, jpl, zs0ice, zs0sn, zs0a, zs0c0 , zs0sm , zs0oi ) CALL wrk_dealloc( jpi, jpj, nlay_i+1, jpl, zs0e ) CALL wrk_dealloc( jpi, jpj, jpl, zviold, zvsold, zaiold, zhimax ) ! clem ! IF( nn_timing == 1 ) CALL timing_stop('limtrp') END SUBROUTINE lim_trp #else !!---------------------------------------------------------------------- !! Default option Empty Module No sea-ice model !!---------------------------------------------------------------------- CONTAINS SUBROUTINE lim_trp ! Empty routine END SUBROUTINE lim_trp #endif !!====================================================================== END MODULE limtrp