MODULE domvvl !!====================================================================== !! *** MODULE domvvl *** !! Ocean : !!====================================================================== !! History : 2.0 ! 2006-06 (B. Levier, L. Marie) original code !! 3.1 ! 2009-02 (G. Madec, M. Leclair, R. Benshila) pure z* coordinate !! 3.3 ! 2011-10 (M. Leclair) totally rewrote domvvl: !! vvl option includes z_star and z_tilde coordinates #if defined key_vvl !!---------------------------------------------------------------------- !! 'key_vvl' variable volume !!---------------------------------------------------------------------- !!---------------------------------------------------------------------- !! dom_vvl_init : define initial vertical scale factors, depths and column thickness !! dom_vvl_nxt : Compute next vertical scale factors !! dom_vvl_swp : Swap vertical scale factors and update the vertical grid !! dom_vvl_interpol : Interpolate vertical scale factors from one grid point to another !! dom_vvl_rst : read/write restart file !! dom_vvl_ctl : Check the vvl options !!---------------------------------------------------------------------- !! * Modules used USE oce ! ocean dynamics and tracers USE dom_oce ! ocean space and time domain USE sbc_oce ! ocean surface boundary condition USE in_out_manager ! I/O manager USE iom ! I/O manager library USE restart ! ocean restart USE lib_mpp ! distributed memory computing library USE lbclnk ! ocean lateral boundary conditions (or mpp link) IMPLICIT NONE PRIVATE !! * Routine accessibility PUBLIC dom_vvl_init ! called by domain.F90 PUBLIC dom_vvl_sf_nxt ! called by step.F90 PUBLIC dom_vvl_sf_swp ! called by step.F90 PUBLIC dom_vvl_interpol ! called by dynnxt.F90 !!* Namelist nam_vvl LOGICAL , PUBLIC :: ln_vvl_zstar = .TRUE. ! zstar vertical coordinate LOGICAL , PUBLIC :: ln_vvl_ztilde = .FALSE. ! ztilde vertical coordinate LOGICAL , PUBLIC :: ln_vvl_layer = .FALSE. ! level vertical coordinate LOGICAL , PUBLIC :: ln_vvl_kepe = .FALSE. ! kinetic/potential energy transfer ! ! conservation: not used yet REAL(wp), PUBLIC :: ahe3 = 0.e0 ! thickness diffusion coefficient !! * Module variables INTEGER :: nvvl ! choice of vertical coordinate REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: un_td, vn_td ! thickness diffusion transport REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: hdiv_lf ! low frequency part of hz divergence REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: e3t_t_b, e3t_t_n, e3t_t_a ! baroclinic scale factors REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:) :: frq_restore_e3t ! retoring period for scale factors REAL(wp) , ALLOCATABLE, SAVE, DIMENSION(:,:) :: frq_restore_hdv ! retoring period for low freq. divergence !! * Substitutions # include "domzgr_substitute.h90" # include "vectopt_loop_substitute.h90" !!---------------------------------------------------------------------- !! NEMO/OPA 3.3 , NEMO-Consortium (2010) !! $Id$ !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) !!---------------------------------------------------------------------- CONTAINS INTEGER FUNCTION dom_vvl_alloc() !!---------------------------------------------------------------------- !! *** FUNCTION dom_vvl_alloc *** !!---------------------------------------------------------------------- IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN ALLOCATE( e3t_t_n(jpi,jpj,jpk) , e3t_t_a(jpi,jpj,jpk) , e3t_t_b(jpi,jpj,jpk) , & & un_td (jpi,jpj,jpk) , vn_td (jpi,jpj,jpk) , STAT = dom_vvl_alloc ) IF( lk_mpp ) CALL mpp_sum ( dom_vvl_alloc ) IF( dom_vvl_alloc /= 0 ) CALL ctl_warn('dom_vvl_alloc: failed to allocate arrays') ENDIF IF( ln_vvl_ztilde ) THEN ALLOCATE( frq_restore_e3t(jpi,jpj) , frq_restore_hdv(jpi,jpj) , hdiv_lf(jpi,jpj,jpk) , STAT= dom_vvl_alloc ) IF( lk_mpp ) CALL mpp_sum ( dom_vvl_alloc ) IF( dom_vvl_alloc /= 0 ) CALL ctl_warn('dom_vvl_alloc: failed to allocate arrays') ENDIF END FUNCTION dom_vvl_alloc SUBROUTINE dom_vvl_init !!---------------------------------------------------------------------- !! *** ROUTINE dom_vvl_init *** !! !! ** Purpose : Initialization of all scale factors, depths !! and water column heights !! !! ** Method : - use restart file and/or initialize !! - interpolate scale factors !! !! ** Action : - fse3t_(n/b) and e3t_t_(n/b) !! - Regrid: fse3(u/v)_n !! fse3(u/v)_b !! fse3w_n !! fse3(u/v)w_b !! fse3(u/v)w_n !! fsdept_n, fsdepw_n and fsde3w_n !! h(u/v) and h(u/v)r !! - ht_0 and ht1_0 !! !! Reference : Leclair, M., and G. Madec, 2011, Ocean Modelling. !!---------------------------------------------------------------------- !! * Local declarations INTEGER :: jk !!---------------------------------------------------------------------- IF(lwp) WRITE(numout,*) IF(lwp) WRITE(numout,*) 'dom_vvl_init : Variable volume activated' IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~' #if defined key_zco CALL ctl_stop( 'dom_vvl_init : options key_zco/key_dynspg_rl are incompatible with variable volume option key_vvl') #endif ! choose vertical coordinate (z_star, z_tilde or layer) ! ========================== CALL dom_vvl_ctl ! Allocate module arrays ! ====================== IF( dom_vvl_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'dom_vvl_init : unable to allocate arrays' ) ! read or initialize e3t_t_(b/n) and fse3t_(b/n) ! ============================================== CALL dom_vvl_rst( nit000, 'READ' ) ! Reconstruction of all vertical scale factors at now and before time steps ! ========================================================================= ! Horizontal scale factor interpolations ! -------------------------------------- CALL dom_vvl_interpol( fse3t_b(:,:,:), fse3u_b(:,:,:), 'U' ) CALL dom_vvl_interpol( fse3t_b(:,:,:), fse3v_b(:,:,:), 'V' ) CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3u_n(:,:,:), 'U' ) CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3v_n(:,:,:), 'V' ) CALL dom_vvl_interpol( fse3u_n(:,:,:), fse3f_n(:,:,:), 'F' ) ! Vertical scale factor interpolations ! ------------------------------------ CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3w_n (:,:,:), 'W' ) CALL dom_vvl_interpol( fse3u_n(:,:,:), fse3uw_n(:,:,:), 'UW' ) CALL dom_vvl_interpol( fse3v_n(:,:,:), fse3vw_n(:,:,:), 'VW' ) CALL dom_vvl_interpol( fse3u_b(:,:,:), fse3uw_b(:,:,:), 'UW' ) CALL dom_vvl_interpol( fse3v_b(:,:,:), fse3vw_b(:,:,:), 'VW' ) ! t- and w- points depth ! ---------------------- fsdept_n(:,:,1) = 0.5 * fse3w_n(:,:,1) fsdepw_n(:,:,1) = 0.e0 fsde3w_n(:,:,1) = fsdept_n(:,:,1) - sshn(:,:) DO jk = 2, jpk fsdept_n(:,:,jk) = fsdept_n(:,:,jk-1) + fse3w_n(:,:,jk) fsdepw_n(:,:,jk) = fsdepw_n(:,:,jk-1) + fse3t_n(:,:,jk-1) fsde3w_n(:,:,jk) = fsdept_n(:,:,jk ) - sshn (:,:) END DO ! Reference water column height at t-, u- and v- point ! ---------------------------------------------------- ht_0(:,:) = 0.e0 hu_0(:,:) = 0.e0 hv_0(:,:) = 0.e0 DO jk = 1, jpk ht_0(:,:) = ht_0(:,:) + fse3t_0(:,:,jk) * tmask(:,:,jk) hu_0(:,:) = hu_0(:,:) + fse3u_0(:,:,jk) * umask(:,:,jk) hv_0(:,:) = hv_0(:,:) + fse3v_0(:,:,jk) * vmask(:,:,jk) END DO ! Restoring frequencies for z_tilde coordinate ! ============================================ IF( ln_vvl_ztilde ) THEN ! - ML - In the future, this should be tunable by the user ! DO jj = 1, jpj ! DO ji = 1, jpi ! frq_restore_hdv(ji,jj) = 2.e0 * rpi / 86400.e0 / 5.e0 & ! & * MAX( SIN( ABS( gphit(ji,jj) ) / .5e0, 1.e0 / 6.e0) ) ! END DO ! END DO ! frq_restore_e3t(:,:) = ( 1.e0 / 6.e0 ) * frq_restore_hdv(:,:) frq_restore_e3t(:,:) = 2.e0 * rpi / ( 30.e0 * 86400.e0 ) frq_restore_hdv(:,:) = 2.e0 * rpi / ( 5.e0 * 86400.e0 ) ! frq_restore_hdv(:,:) = 2.e0 * rpi / ( 2.e0 * 86400.e0 ) ENDIF END SUBROUTINE dom_vvl_init SUBROUTINE dom_vvl_sf_nxt( kt ) !!---------------------------------------------------------------------- !! *** ROUTINE dom_vvl_sf_nxt *** !! !! ** Purpose : - compute the after scale factors used in tra_zdf, dynnxt, !! tranxt and dynspg routines !! !! ** Method : - z_tilde_case: after scale factor increment computed with !! high frequency part of horizontal divergence + retsoring to !! towards the background grid + thickness difusion. !! Then repartition of ssh INCREMENT proportionnaly !! to the "baroclinic" level thickness. !! - z_star case: Repartition of ssh proportionnaly to the level thickness. !! !! ** Action : - hdiv_lf: restoring towards full baroclinic divergence in z_tilde case !! - e3t_t_a: after increment of vertical scale factor !! in z_tilde case !! - fse3(t/u/v)_a !! !! Reference : Leclair, M., and G. Madec, 2011, Ocean Modelling. !!---------------------------------------------------------------------- USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released USE oce , ONLY: ze3t => ta ! ua used as workspace USE wrk_nemo, ONLY: zht => wrk_2d_1 ! 2D REAL workspace USE wrk_nemo, ONLY: z_scale => wrk_2d_2 ! 2D REAL workspace USE wrk_nemo, ONLY: zwu => wrk_2d_3 ! 2D REAL workspace USE wrk_nemo, ONLY: zwv => wrk_2d_4 ! 2D REAL workspace USE wrk_nemo, ONLY: zhdiv => wrk_2d_5 ! 2D REAL workspace !! * Arguments INTEGER, INTENT( in ) :: kt ! time step !! * Local declarations INTEGER :: ji, jj, jk ! dummy loop indices INTEGER , DIMENSION(3) :: ijk_max, ijk_min ! temporary integers REAL(wp) :: z2dt ! temporary scalars REAL(wp) :: z_def_max ! temporary scalar REAL(wp) :: z_tmin, z_tmax ! temporary scalars LOGICAL :: ln_debug = .FALSE. ! local logical for debug prints !!---------------------------------------------------------------------- IF( wrk_in_use(2, 1, 2, 3, 4, 5) ) THEN CALL ctl_stop('dom_vvl_sf_nxt: requested workspace arrays unavailable') ; RETURN END IF IF(kt == nit000) THEN IF(lwp) WRITE(numout,*) IF(lwp) WRITE(numout,*) 'dom_vvl_sf_nxt : compute after scale factors' IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~' ENDIF ! ******************************* ! ! After acale factors at t-points ! ! ******************************* ! ! ! ----------------- ! IF( ln_vvl_zstar ) THEN ! z_star coordinate ! ! ! ----------------- ! z_scale(:,:) = ( ssha(:,:) - sshb(:,:) ) * tmask(:,:,1) / ( ht_0(:,:) + sshn(:,:) + 1. - tmask(:,:,1) ) DO jk = 1, jpkm1 fse3t_a(:,:,jk) = fse3t_b(:,:,jk) + fse3t_n(:,:,jk) * z_scale(:,:) END DO ! ! --------------------------- ! ELSE ! z_tilde or layer coordinate ! ! ! --------------------------- ! ! I - Low frequency horizontal divergence ! ======================================= ! 1 - barotropic divergence ! ------------------------- zhdiv(:,:) = 0. zht(:,:) = 0. DO jk = 1, jpkm1 zhdiv(:,:) = zhdiv(:,:) + hdivn(:,:,jk) * fse3t_n(:,:,jk) zht (:,:) = zht (:,:) + fse3t_n(:,:,jk) * tmask(:,:,jk) END DO zhdiv(:,:) = zhdiv(:,:) / ( zht(:,:) + 1. - tmask(:,:,1) ) ! 2 - restoring equation (z-tilde case only) ! ---------------------- IF( ln_vvl_ztilde ) THEN IF( kt .GT. nit000 ) THEN DO jk = 1, jpkm1 hdiv_lf(:,:,jk) = hdiv_lf(:,:,jk) - rdt * frq_restore_hdv(:,:) & & * ( hdiv_lf(:,:,jk) - fse3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) ) END DO ENDIF END IF ! II - after z_tilde increments of vertical scale factors ! ======================================================= e3t_t_a(:,:,:) = 0.e0 ! e3t_t_a used to store tendency terms ! 1 - High frequency divergence term ! ---------------------------------- IF( ln_vvl_ztilde ) THEN ! z_tilde case DO jk = 1, jpkm1 e3t_t_a(:,:,jk) = e3t_t_a(:,:,jk) - ( fse3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) - hdiv_lf(:,:,jk) ) END DO ! layer case ELSE DO jk = 1, jpkm1 e3t_t_a(:,:,jk) = e3t_t_a(:,:,jk) - fse3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) END DO END IF ! 2 - Restoring term (z-tilde case only) ! ------------------ IF( ln_vvl_ztilde ) THEN DO jk = 1, jpk e3t_t_a(:,:,jk) = e3t_t_a(:,:,jk) - frq_restore_e3t(:,:) * e3t_t_b(:,:,jk) END DO END IF ! 3 - Thickness diffusion term ! ---------------------------- zwu(:,:) = 0.e0 zwv(:,:) = 0.e0 ! a - first derivative: diffusive fluxes DO jk = 1, jpkm1 DO jj = 1, jpjm1 DO ji = 1, fs_jpim1 ! vector opt. un_td(ji,jj,jk) = ahe3 * umask(ji,jj,jk) * e1ur(ji,jj) * ( e3t_t_b(ji,jj,jk) - e3t_t_b(ji+1,jj ,jk) ) vn_td(ji,jj,jk) = ahe3 * vmask(ji,jj,jk) * e2vr(ji,jj) * ( e3t_t_b(ji,jj,jk) - e3t_t_b(ji ,jj+1,jk) ) zwu(ji,jj) = zwu(ji,jj) + un_td(ji,jj,jk) zwv(ji,jj) = zwv(ji,jj) + vn_td(ji,jj,jk) END DO END DO END DO ! b - correction for last oceanic u-v points DO jj = 1, jpj DO ji = 1, jpi un_td(ji,jj,mbku(ji,jj)) = un_td(ji,jj,mbku(ji,jj)) - zwu(ji,jj) vn_td(ji,jj,mbkv(ji,jj)) = vn_td(ji,jj,mbkv(ji,jj)) - zwv(ji,jj) END DO END DO ! c - second derivative: divergence of diffusive fluxes DO jk = 1, jpkm1 DO jj = 2, jpjm1 DO ji = fs_2, fs_jpim1 ! vector opt. e3t_t_a(ji,jj,jk) = e3t_t_a(ji,jj,jk) + ( un_td(ji-1,jj ,jk) - un_td(ji,jj,jk) & & + vn_td(ji ,jj-1,jk) - vn_td(ji,jj,jk) ) & & * e12t_1(ji,jj) END DO END DO END DO ! d - thickness diffusion equivalent transport: boundary conditions ! (stored for tracer advction and continuity equation) CALL lbc_lnk( un_td , 'U' , -1.) CALL lbc_lnk( vn_td , 'V' , -1.) ! 4 - Time stepping of baroclinic scale factors ! --------------------------------------------- ! Leapfrog time stepping ! ~~~~~~~~~~~~~~~~~~~~~~ IF( neuler == 0 .AND. kt == nit000 ) THEN z2dt = rdt ELSE z2dt = 2.e0 * rdt ENDIF CALL lbc_lnk( e3t_t_a(:,:,:), 'T', 1. ) e3t_t_a(:,:,:) = e3t_t_b(:,:,:) + z2dt * tmask(:,:,:) * e3t_t_a(:,:,:) fse3t_a(:,:,:) = fse3t_0(:,:,:) + e3t_t_a(:,:,:) ! Maximum deformation control ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~ ! - ML - Should perhaps be put in the namelist z_def_max = 0.9e0 ze3t(:,:,jpk) = 0.e0 DO jk = 1, jpkm1 ze3t(:,:,jk) = e3t_t_a(:,:,jk) / fse3t_0(:,:,jk) * tmask(:,:,jk) * tmask_i(:,:) END DO z_tmax = MAXVAL( ze3t(:,:,:) ) z_tmin = MINVAL( ze3t(:,:,:) ) ! - ML - test: for the moment, stop simulation for too large e3_t variations IF( ( z_tmax .GT. z_def_max ) .OR. ( z_tmin .LT. - z_def_max ) ) THEN ijk_max = MAXLOC( ze3t(:,:,:) ) ijk_min = MINLOC( ze3t(:,:,:) ) WRITE(numout, *) 'MAX( e3t_t_a(:,:,:) / fse3t_0(:,:,:) ) =', z_tmax WRITE(numout, *) 'at i, j, k=', ijk_max WRITE(numout, *) 'MIN( e3t_t_a(:,:,:) / fse3t_0(:,:,:) ) =', z_tmin WRITE(numout, *) 'at i, j, k=', ijk_min CALL ctl_stop('MAX( ABS( e3t_t_a(:,:,:) ) / fse3t_0(:,:,:) ) too high') ENDIF ! - ML - end test ! - ML - This will cause a baroclinicity error if the ctl_stop above is not used e3t_t_a(:,:,:) = MIN( e3t_t_a(:,:,:), ( 1.e0 + z_def_max ) * fse3t_0(:,:,:) ) e3t_t_a(:,:,:) = MAX( e3t_t_a(:,:,:), ( 1.e0 - z_def_max ) * fse3t_0(:,:,:) ) ! III - Barotropic repartition of the sea surface height over the baroclinic profile ! ================================================================================== ! add e3t(n-1) "star" Asselin-filtered DO jk = 1, jpkm1 fse3t_a(:,:,jk) = fse3t_a(:,:,jk) + fse3t_b(:,:,jk) - fse3t_0(:,:,jk) - e3t_t_b(:,:,jk) END DO ! add ( ssh increment + "baroclinicity error" ) proportionnaly to e3t(n) ! - ML - baroclinicity error should be better treated in the future ! i.e. locally and not spread over the water column. zht(:,:) = 0. DO jk = 1, jpkm1 zht(:,:) = zht(:,:) + e3t_t_a(:,:,jk) * tmask(:,:,jk) END DO z_scale(:,:) = ( ssha(:,:) - sshb(:,:) - zht(:,:) ) / ( ht_0(:,:) + sshn(:,:) + 1. - tmask(:,:,1) ) DO jk = 1, jpkm1 fse3t_a(:,:,jk) = fse3t_a(:,:,jk) + fse3t_n(:,:,jk) * z_scale(:,:) * tmask(:,:,jk) END DO ENDIF IF( ln_debug ) THEN ! - ML - test: control prints for debuging IF ( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN WRITE(numout, *) 'kt =', kt WRITE(numout, *) 'MAXVAL(abs(ht_0-SUM(e3t_t_a))) =', & & MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) - zht(:,:) ) ) END IF zht(:,:) = 0.e0 DO jk = 1, jpkm1 zht(:,:) = zht(:,:) + fse3t_n(:,:,jk) * tmask(:,:,jk) END DO WRITE(numout, *) 'MAXVAL(abs(ht_0+sshn-SUM(fse3t_n))) =', & & MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + sshn(:,:) - zht(:,:) ) ) zht(:,:) = 0.e0 DO jk = 1, jpkm1 zht(:,:) = zht(:,:) + fse3t_a(:,:,jk) * tmask(:,:,jk) END DO WRITE(numout, *) 'MAXVAL(abs(ht_0+sshn-SUM(fse3t_a))) =', & & MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + ssha(:,:) - zht(:,:) ) ) END IF ! *********************************** ! ! After scale factors at u- v- points ! ! *********************************** ! CALL dom_vvl_interpol( fse3t_a(:,:,:), fse3u_a(:,:,:), 'U' ) CALL dom_vvl_interpol( fse3t_a(:,:,:), fse3v_a(:,:,:), 'V' ) IF( wrk_not_released(2, 1, 2, 3, 4, 5) ) THEN CALL ctl_stop( 'dom_vvl_sf_nxt: failed to release workspace arrays' ) ENDIF END SUBROUTINE dom_vvl_sf_nxt SUBROUTINE dom_vvl_sf_swp( kt ) !!---------------------------------------------------------------------- !! *** ROUTINE dom_vvl_sf_swp *** !! !! ** Purpose : compute time filter and swap of scale factors !! compute all depths and related variables for next time step !! write outputs and restart file !! !! ** Method : - swap of e3t with trick for volume/tracer conservation !! - reconstruct scale factor at other grid points (interpolate) !! - recompute depths and water height fields !! !! ** Action : - fse3t_(b/n), e3t_t_(b/n) and fse3(u/v)_n ready for next time step !! - Recompute: !! fse3(u/v)_b !! fse3w_n !! fse3(u/v)w_b !! fse3(u/v)w_n !! fsdept_n, fsdepw_n and fsde3w_n !! h(u/v) and h(u/v)r !! !! Reference : Leclair, M., and G. Madec, 2009, Ocean Modelling. !! Leclair, M., and G. Madec, 2011, Ocean Modelling. !!---------------------------------------------------------------------- USE oce, ONLY: z_e3t_def => ta ! square of baroclinic scale factors deformation (in %) !! * Arguments INTEGER, INTENT( in ) :: kt ! time step !! * Local declarations INTEGER :: jk ! dummy loop indices !!---------------------------------------------------------------------- IF( kt == nit000 ) THEN IF(lwp) WRITE(numout,*) IF(lwp) WRITE(numout,*) 'dom_vvl_sf_swp : - time filter and swap of scale factors' IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~ - interpolate scale factors and compute depths for next time step' ENDIF ! Time filter and swap of scale factors ! ===================================== ! - ML - fse3(t/u/v)_b are allready computed in dynnxt. IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN IF( neuler == 0 .AND. kt == nit000 ) THEN e3t_t_b(:,:,:) = e3t_t_n(:,:,:) ELSE e3t_t_b(:,:,:) = e3t_t_n(:,:,:) + atfp * (e3t_t_b(:,:,:) - 2.e0 * e3t_t_n(:,:,:) + e3t_t_a(:,:,:) ) ENDIF e3t_t_n(:,:,:) = e3t_t_a(:,:,:) ENDIF fse3t_n(:,:,:) = fse3t_a(:,:,:) fse3u_n(:,:,:) = fse3u_a(:,:,:) fse3v_n(:,:,:) = fse3v_a(:,:,:) ! Compute all missing vertical scale factor and depths ! ==================================================== ! Horizontal scale factor interpolations ! -------------------------------------- ! - ML - fse3u_b and fse3v_b are allready computed in dynnxt CALL dom_vvl_interpol( fse3u_n(:,:,:), fse3f_n (:,:,:), 'F' ) ! Vertical scale factor interpolations ! ------------------------------------ CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3w_n (:,:,:), 'W' ) CALL dom_vvl_interpol( fse3u_n(:,:,:), fse3uw_n(:,:,:), 'UW' ) CALL dom_vvl_interpol( fse3v_n(:,:,:), fse3vw_n(:,:,:), 'VW' ) CALL dom_vvl_interpol( fse3u_b(:,:,:), fse3uw_b(:,:,:), 'UW' ) CALL dom_vvl_interpol( fse3v_b(:,:,:), fse3vw_b(:,:,:), 'VW' ) ! t- and w- points depth ! ---------------------- fsdept_n(:,:,1) = 0.5 * fse3w_n(:,:,1) fsdepw_n(:,:,1) = 0.e0 fsde3w_n(:,:,1) = fsdept_n(:,:,1) - sshn(:,:) DO jk = 2, jpk fsdept_n(:,:,jk) = fsdept_n(:,:,jk-1) + fse3w_n(:,:,jk) fsdepw_n(:,:,jk) = fsdepw_n(:,:,jk-1) + fse3t_n(:,:,jk-1) fsde3w_n(:,:,jk) = fsdept_n(:,:,jk ) - sshn (:,:) END DO ! Local depth and Inverse of the local depth of the water column at u- and v- points ! ---------------------------------------------------------------------------------- hu(:,:) = 0. hv(:,:) = 0. DO jk = 1, jpk hu(:,:) = hu(:,:) + fse3u_n(:,:,jk) * umask(:,:,jk) hv(:,:) = hv(:,:) + fse3v_n(:,:,jk) * vmask(:,:,jk) END DO ! Inverse of the local depth hur(:,:) = umask(:,:,1) / ( hu(:,:) + 1. - umask(:,:,1) ) hvr(:,:) = vmask(:,:,1) / ( hv(:,:) + 1. - umask(:,:,1) ) ! Write outputs ! ============= ! - ML - add output variables in xml file for all configurations z_e3t_def(:,:,:) = ( ( fse3t_n(:,:,:) - fse3t_0(:,:,:) ) / fse3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 CALL iom_put( "e3t_n" , fse3t_n (:,:,:) ) CALL iom_put( "dept_n" , fsdept_n (:,:,:) ) CALL iom_put( "e3tdef" , z_e3t_def(:,:,:) ) ! write restart file ! ================== IF( lrst_oce ) CALL dom_vvl_rst( kt, 'WRITE' ) END SUBROUTINE dom_vvl_sf_swp SUBROUTINE dom_vvl_interpol( pe3_in, pe3_out, pout ) !!--------------------------------------------------------------------- !! *** ROUTINE dom_vvl__interpol *** !! !! ** Purpose : interpolate scale factors from one grid point to another !! !! ** Method : e3_out = e3_0 + interpolation(e3_in - e3_0) !! - horizontal interpolation: grid cell surface averaging !! - vertical interpolation: simple averaging !!---------------------------------------------------------------------- !! * Arguments REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( in ) :: pe3_in ! input e3 to be interpolated REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( inout ) :: pe3_out ! output interpolated e3 CHARACTER(len=1), INTENT( in ) :: pout ! grid point of out scale factors ! ! = 'U', 'V', 'W, 'F', 'UW' or 'VW' !! * Local declarations INTEGER :: ji, jj, jk ! dummy loop indices !!---------------------------------------------------------------------- SELECT CASE ( pout ) ! ! ------------------------------------- ! CASE( 'U' ) ! interpolation from T-point to U-point ! ! ! ------------------------------------- ! ! horizontal surface weighted interpolation DO jk = 1, jpkm1 DO jj = 2, jpjm1 DO ji = 1, fs_jpim1 ! vector opt. pe3_out(ji,jj,jk) = umask(ji,jj,jk) * e12u_1(ji,jj) & & * ( e12t(ji ,jj) * ( pe3_in(ji ,jj,jk) - fse3t_0(ji ,jj,jk) ) & & + e12t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - fse3t_0(ji+1,jj,jk) ) ) END DO END DO END DO ! boundary conditions CALL lbc_lnk( pe3_out(:,:,:), 'U', 1. ) pe3_out(:,:,:) = pe3_out(:,:,:) + fse3u_0(:,:,:) ! ! ------------------------------------- ! CASE( 'V' ) ! interpolation from T-point to V-point ! ! ! ------------------------------------- ! ! horizontal surface weighted interpolation DO jk = 1, jpkm1 DO jj = 1, jpjm1 DO ji = fs_2, fs_jpim1 ! vector opt. pe3_out(ji,jj,jk) = umask(ji,jj,jk) * e12v_1(ji,jj) & & * ( e12t(ji,jj ) * ( pe3_in(ji,jj ,jk) - fse3t_0(ji,jj ,jk) ) & & + e12t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - fse3t_0(ji,jj+1,jk) ) ) END DO END DO END DO ! boundary conditions CALL lbc_lnk( pe3_out(:,:,:), 'V', 1. ) pe3_out(:,:,:) = pe3_out(:,:,:) + fse3v_0(:,:,:) ! ! ------------------------------------- ! CASE( 'F' ) ! interpolation from U-point to F-point ! ! ! ------------------------------------- ! ! horizontal surface weighted interpolation DO jk = 1, jpkm1 DO jj = 1, jpjm1 DO ji = 1, fs_jpim1 ! vector opt. pe3_out(ji,jj,jk) = umask(ji,jj,jk) * umask(ji,jj+1,jk) * e12f_1(ji,jj) & & * ( e12u(ji,jj ) * ( pe3_in(ji,jj ,jk) - fse3u_0(ji,jj ,jk) ) & & + e12u(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - fse3u_0(ji,jj+1,jk) ) ) END DO END DO END DO ! boundary conditions CALL lbc_lnk( pe3_out(:,:,:), 'F', 1. ) pe3_out(:,:,:) = pe3_out(:,:,:) + fse3f_0(:,:,:) ! ! ------------------------------------- ! CASE( 'W' ) ! interpolation from T-point to W-point ! ! ! ------------------------------------- ! ! vertical simple interpolation pe3_out(:,:,1) = fse3w_0(:,:,1) + pe3_in(:,:,1) - fse3t_0(:,:,1) ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without undirect adressing DO jk = 2, jpk pe3_out(:,:,jk) = fse3w_0(:,:,jk) + ( 1. - 0.5 * tmask(:,:,jk) ) * ( pe3_in(:,:,jk-1) - fse3t_0(:,:,jk-1) ) & & + ( 0.5 * tmask(:,:,jk) ) * ( pe3_in(:,:,jk ) - fse3t_0(:,:,jk ) ) END DO ! ! -------------------------------------- ! CASE( 'UW' ) ! interpolation from U-point to UW-point ! ! ! -------------------------------------- ! ! vertical simple interpolation pe3_out(:,:,1) = fse3uw_0(:,:,1) + pe3_in(:,:,1) - fse3u_0(:,:,1) ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without undirect adressing DO jk = 2, jpk pe3_out(:,:,jk) = fse3uw_0(:,:,jk) + ( 1. - 0.5 * umask(:,:,jk) ) * ( pe3_in(:,:,jk-1) - fse3u_0(:,:,jk-1) ) & & + ( 0.5 * umask(:,:,jk) ) * ( pe3_in(:,:,jk ) - fse3u_0(:,:,jk ) ) END DO ! ! -------------------------------------- ! CASE( 'VW' ) ! interpolation from V-point to VW-point ! ! ! -------------------------------------- ! ! vertical simple interpolation pe3_out(:,:,1) = fse3vw_0(:,:,1) + pe3_in(:,:,1) - fse3v_0(:,:,1) ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without undirect adressing DO jk = 2, jpk pe3_out(:,:,jk) = fse3vw_0(:,:,jk) + ( 1. - 0.5 * vmask(:,:,jk) ) * ( pe3_in(:,:,jk-1) - fse3v_0(:,:,jk-1) ) & & + ( 0.5 * vmask(:,:,jk) ) * ( pe3_in(:,:,jk ) - fse3v_0(:,:,jk ) ) END DO END SELECT END SUBROUTINE dom_vvl_interpol SUBROUTINE dom_vvl_rst( kt, cdrw ) !!--------------------------------------------------------------------- !! *** ROUTINE dom_vvl_rst *** !! !! ** Purpose : Read or write VVL file in restart file !! !! ** Method : use of IOM library !! if the restart does not contain vertical scale factors, !! they are set to the _0 values !! if the restart does not contain vertical scale factors increments (z_tilde), !! they are set to 0. !!---------------------------------------------------------------------- !! * Arguments INTEGER , INTENT(in) :: kt ! ocean time-step CHARACTER(len=*), INTENT(in) :: cdrw ! "READ"/"WRITE" flag !! * Local declarations INTEGER :: id1, id2, id3, id4, id5 ! local integers !!---------------------------------------------------------------------- ! IF( TRIM(cdrw) == 'READ' ) THEN ! Read/initialise ! ! =============== IF( ln_rstart ) THEN !* Read the restart file id1 = iom_varid( numror, 'fse3t_b', ldstop = .FALSE. ) id2 = iom_varid( numror, 'fse3t_n', ldstop = .FALSE. ) id3 = iom_varid( numror, 'e3t_t_b', ldstop = .FALSE. ) id4 = iom_varid( numror, 'e3t_t_n', ldstop = .FALSE. ) id5 = iom_varid( numror, 'hdif_lf', ldstop = .FALSE. ) ! ! --------- ! ! ! all cases ! ! ! --------- ! IF( MIN( id1, id2 ) > 0 ) THEN ! all required arrays exist CALL iom_get( numror, jpdom_autoglo, 'fse3t_b', fse3t_b(:,:,:) ) CALL iom_get( numror, jpdom_autoglo, 'fse3t_n', fse3t_n(:,:,:) ) IF( neuler == 0 ) THEN fse3t_b(:,:,:) = fse3t_n(:,:,:) ENDIF ELSE ! one at least array is missing CALL ctl_stop( 'dom_vvl_rst: vvl cannot restart from a non vvl run' ) ENDIF ! ! ----------- ! IF( ln_vvl_zstar ) THEN ! z_star case ! ! ! ----------- ! IF( MIN( id3, id4 ) > 0 ) THEN CALL ctl_stop( 'dom_vvl_rst: z_star cannot restart from a z_tilde or layer run' ) ENDIF ! ! ----------------------- ! ELSE ! z_tilde and layer cases ! ! ! ----------------------- ! IF( MIN( id3, id4 ) > 0 ) THEN ! all required arrays exist CALL iom_get( numror, jpdom_autoglo, 'e3t_t_b', e3t_t_b(:,:,:) ) CALL iom_get( numror, jpdom_autoglo, 'e3t_t_n', e3t_t_n(:,:,:) ) ELSE ! one at least array is missing e3t_t_b(:,:,:) = 0.e0 e3t_t_n(:,:,:) = 0.e0 ENDIF ! ! ------------ ! IF( ln_vvl_ztilde ) THEN ! z_tilde case ! ! ! ------------ ! IF( id5 > 0 ) THEN ! required array exists CALL iom_get( numror, jpdom_autoglo, 'hdiv_lf', hdiv_lf(:,:,:) ) ELSE ! array is missing hdiv_lf(:,:,:) = 0.e0 ENDIF ENDIF ENDIF ! ELSE !* Initialize at "rest" fse3t_b(:,:,:) = fse3t_0(:,:,:) fse3t_n(:,:,:) = fse3t_0(:,:,:) e3t_t_b(:,:,:) = 0.e0 e3t_t_n(:,:,:) = 0.e0 hdiv_lf(:,:,:) = 0.e0 ENDIF ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN ! Create restart file ! ! =================== IF(lwp) WRITE(numout,*) '---- dom_vvl_rst ----' ! ! --------- ! ! ! all cases ! ! ! --------- ! CALL iom_rstput( kt, nitrst, numrow, 'fse3t_b', fse3t_b(:,:,:) ) CALL iom_rstput( kt, nitrst, numrow, 'fse3t_n', fse3t_n(:,:,:) ) ! ! ----------------------- ! IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN ! z_tilde and layer cases ! ! ! ----------------------- ! CALL iom_rstput( kt, nitrst, numrow, 'e3t_t_b', e3t_t_b(:,:,:) ) CALL iom_rstput( kt, nitrst, numrow, 'e3t_t_n', e3t_t_n(:,:,:) ) END IF ! ! -------------! IF( ln_vvl_ztilde ) THEN ! z_tilde case ! ! ! ------------ ! CALL iom_rstput( kt, nitrst, numrow, 'hdiv_lf', hdiv_lf(:,:,:) ) ENDIF ENDIF END SUBROUTINE dom_vvl_rst SUBROUTINE dom_vvl_ctl !!--------------------------------------------------------------------- !! *** ROUTINE dom_vvl_ctl *** !! !! ** Purpose : Control the consistency between namelist options for !! vertical coordinate and set nvvl !!---------------------------------------------------------------------- INTEGER :: ioptio NAMELIST/nam_vvl/ ln_vvl_zstar, ln_vvl_ztilde, ln_vvl_layer, ahe3! , ln_vvl_kepe !!---------------------------------------------------------------------- REWIND ( numnam ) ! Read Namelist nam_vvl : vertical coordinate READ ( numnam, nam_vvl ) IF(lwp) THEN ! Namelist print WRITE(numout,*) WRITE(numout,*) 'dom_vvl_ctl : choice/control of the variable vertical coordinate' WRITE(numout,*) '~~~~~~~~~~~' WRITE(numout,*) ' Namelist nam_vvl : chose a vertical coordinate' WRITE(numout,*) ' zstar ln_vvl_zstar = ', ln_vvl_zstar WRITE(numout,*) ' ztilde ln_vvl_ztilde = ', ln_vvl_ztilde WRITE(numout,*) ' layer ln_vvl_layer = ', ln_vvl_layer ! WRITE(numout,*) ' Namelist nam_vvl : chose kinetic-to-potential energy conservation' ! WRITE(numout,*) ' ln_vvl_kepe = ', ln_vvl_kepe WRITE(numout,*) ' Namelist nam_vvl : thickness diffusion coefficient' WRITE(numout,*) ' ahe3 = ', ahe3 ENDIF ioptio = 0 ! Parameter control IF( ln_vvl_zstar ) ioptio = ioptio + 1 IF( ln_vvl_ztilde ) ioptio = ioptio + 1 IF( ln_vvl_layer ) ioptio = ioptio + 1 IF( ioptio /= 1 ) CALL ctl_stop( 'Choose ONE vertical coordinate in namelist nam_vvl' ) IF( ln_vvl_zstar ) nvvl = 1 IF( ln_vvl_ztilde ) nvvl = 2 IF( ln_vvl_layer ) nvvl = 3 IF(lwp) THEN ! Print the choice WRITE(numout,*) IF( nvvl == 1 ) WRITE(numout,*) ' zstar vertical coordinate is used' IF( nvvl == 2 ) WRITE(numout,*) ' ztilde vertical coordinate is used' IF( nvvl == 3 ) WRITE(numout,*) ' layer vertical coordinate is used' ! - ML - Option not developed yet ! IF( ln_vvl_kepe ) WRITE(numout,*) ' kinetic to potential energy transfer : option used' ! IF( .NOT. ln_vvl_kepe ) WRITE(numout,*) ' kinetic to potential energy transfer : option not used' ENDIF END SUBROUTINE dom_vvl_ctl #else !!---------------------------------------------------------------------- !! Default option : Empty routine !!---------------------------------------------------------------------- SUBROUTINE dom_vvl_init WRITE(*,*) 'dom_vvl_init: You should not have seen this print! error?' END SUBROUTINE dom_vvl_init SUBROUTINE dom_vvl_sf_nxt( kt ) !! * Arguments INTEGER, INTENT( in ) :: kt WRITE(*,*) 'dom_vvl_sf_nxt: You should not have seen this print! error?', kt END SUBROUTINE dom_vvl_sf_nxt SUBROUTINE dom_vvl_sf_swp( kt ) !! * Arguments INTEGER, INTENT( in ) :: kt WRITE(*,*) 'dom_vvl_sf_swp: You should not have seen this print! error?', kt END SUBROUTINE dom_vvl_sf_swp SUBROUTINE dom_vvl_interpol( pe3_in, pe3_out, pout ) REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( in ) :: pe3_in REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( inout ) :: pe3_out CHARACTER(len=1), INTENT( in ) :: pout WRITE(*,*) 'dom_vvl_interpol: You should not have seen this print! error?' END SUBROUTINE dom_vvl_interpol #endif !!====================================================================== END MODULE domvvl