MODULE domqe !!====================================================================== !! *** MODULE domqe *** !! 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 !! 3.6 ! 2014-11 (P. Mathiot) add ice shelf capability !! 4.1 ! 2019-08 (A. Coward, D. Storkey) rename dom_vvl_sf_swp -> dom_vvl_sf_update for new timestepping !! 4.x ! 2020-02 (G. Madec, S. Techene) pure z* (quasi-eulerian) coordinate !!---------------------------------------------------------------------- !!---------------------------------------------------------------------- !! dom_qe_init : define initial vertical scale factors, depths and column thickness !! dom_qe_sf_nxt : Compute next vertical scale factors !! dom_qe_sf_update : Swap vertical scale factors and update the vertical grid !! dom_qe_interpol : Interpolate vertical scale factors from one grid point to another !! dom_qe_r3c : Compute ssh/h_0 ratioat t-, u-, v-, and optionally f-points !! dom_qe_rst : read/write restart file !! dom_qe_ctl : Check the vvl options !!---------------------------------------------------------------------- USE oce ! ocean dynamics and tracers USE phycst ! physical constant USE dom_oce ! ocean space and time domain USE dynadv , ONLY : ln_dynadv_vec USE isf_oce ! iceshelf cavities USE sbc_oce ! ocean surface boundary condition USE wet_dry ! wetting and drying USE usrdef_istate ! user defined initial state (wad only) USE restart ! ocean restart ! USE in_out_manager ! I/O manager USE iom ! I/O manager library USE lib_mpp ! distributed memory computing library USE lbclnk ! ocean lateral boundary conditions (or mpp link) USE timing ! Timing IMPLICIT NONE PRIVATE PUBLIC dom_qe_init ! called by domain.F90 PUBLIC dom_qe_zgr ! called by isfcpl.F90 PUBLIC dom_qe_sf_nxt ! called by steplf.F90 PUBLIC dom_qe_sf_update ! called by steplf.F90 PUBLIC dom_qe_interpol ! called by dynnxt.F90 PUBLIC dom_qe_r3c ! called by steplf.F90 ! !!* Namelist nam_vvl LOGICAL , PUBLIC :: ln_vvl_zstar = .FALSE. ! 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_ztilde_as_zstar = .FALSE. ! ztilde vertical coordinate LOGICAL , PUBLIC :: ln_vvl_zstar_at_eqtor = .FALSE. ! ztilde vertical coordinate LOGICAL , PUBLIC :: ln_vvl_kepe = .FALSE. ! kinetic/potential energy transfer ! ! conservation: not used yet REAL(wp) :: rn_ahe3 ! thickness diffusion coefficient REAL(wp) :: rn_rst_e3t ! ztilde to zstar restoration timescale [days] REAL(wp) :: rn_lf_cutoff ! cutoff frequency for low-pass filter [days] REAL(wp) :: rn_zdef_max ! maximum fractional e3t deformation LOGICAL , PUBLIC :: ln_vvl_dbg = .FALSE. ! debug control prints REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: un_td, vn_td ! thickness diffusion transport !! * Substitutions # include "do_loop_substitute.h90" !!---------------------------------------------------------------------- !! NEMO/OCE 4.0 , NEMO Consortium (2018) !! $Id: domvvl.F90 12377 2020-02-12 14:39:06Z acc $ !! Software governed by the CeCILL license (see ./LICENSE) !!---------------------------------------------------------------------- CONTAINS SUBROUTINE dom_qe_init( Kbb, Kmm, Kaa ) !!---------------------------------------------------------------------- !! *** ROUTINE dom_qe_init *** !! !! ** Purpose : Initialization of all scale factors, depths !! and water column heights !! !! ** Method : - use restart file and/or initialize !! - interpolate scale factors !! !! ** Action : - e3t_(n/b) !! - Regrid: e3[u/v](:,:,:,Kmm) !! e3[u/v](:,:,:,Kmm) !! e3w(:,:,:,Kmm) !! e3[u/v]w_b !! e3[u/v]w_n !! gdept(:,:,:,Kmm), gdepw(:,:,:,Kmm) and gde3w !! - h(t/u/v)_0 !! !! Reference : Leclair, M., and G. Madec, 2011, Ocean Modelling. !!---------------------------------------------------------------------- INTEGER, INTENT(in) :: Kbb, Kmm, Kaa ! IF(lwp) WRITE(numout,*) IF(lwp) WRITE(numout,*) 'dom_qe_init : Variable volume activated' IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~' ! CALL dom_qe_ctl ! choose vertical coordinate (z_star, z_tilde or layer) ! ! ! Read or initialize e3t_(b/n), tilde_e3t_(b/n) and hdiv_lf CALL dom_qe_rst( nit000, Kbb, Kmm, 'READ' ) e3t(:,:,jpk,Kaa) = e3t_0(:,:,jpk) ! last level always inside the sea floor set one for all ! CALL dom_qe_zgr(Kbb, Kmm, Kaa) ! interpolation scale factor, depth and water column ! IF(lwxios) THEN ! define variables in restart file when writing with XIOS CALL iom_set_rstw_var_active('e3t_b') CALL iom_set_rstw_var_active('e3t_n') ENDIF ! END SUBROUTINE dom_qe_init SUBROUTINE dom_qe_zgr(Kbb, Kmm, Kaa) !!---------------------------------------------------------------------- !! *** ROUTINE dom_qe_init *** !! !! ** Purpose : Interpolation of all scale factors, !! depths and water column heights !! !! ** Method : - interpolate scale factors !! !! ** Action : - e3t_(n/b) !! - Regrid: e3(u/v)_n !! e3(u/v)_b !! e3w_n !! e3(u/v)w_b !! e3(u/v)w_n !! gdept_n, gdepw_n and gde3w_n !! - h(t/u/v)_0 !! !! Reference : Leclair, M., and G. Madec, 2011, Ocean Modelling. !!---------------------------------------------------------------------- INTEGER, INTENT(in) :: Kbb, Kmm, Kaa !!---------------------------------------------------------------------- INTEGER :: ji, jj, jk INTEGER :: ii0, ii1, ij0, ij1 REAL(wp):: zcoef !!---------------------------------------------------------------------- ! ! !== Set of all other vertical scale factors ==! (now and before) ! ! Horizontal interpolation of e3t CALL dom_qe_r3c( ssh(:,:,Kbb), r3t(:,:,Kbb), r3u(:,:,Kbb), r3v(:,:,Kbb) ) CALL dom_qe_r3c( ssh(:,:,Kmm), r3t(:,:,Kmm), r3u(:,:,Kmm), r3v(:,:,Kmm), r3f(:,:) ) ! DO jk = 1, jpkm1 ! Horizontal interpolation of e3t e3t(:,:,jk,Kbb) = e3t_0(:,:,jk) * ( 1._wp + r3t(:,:,Kbb) * tmask(:,:,jk) ) ! Kbb time level e3u(:,:,jk,Kbb) = e3u_0(:,:,jk) * ( 1._wp + r3u(:,:,Kbb) * umask(:,:,jk) ) e3v(:,:,jk,Kbb) = e3v_0(:,:,jk) * ( 1._wp + r3v(:,:,Kbb) * vmask(:,:,jk) ) e3t(:,:,jk,Kmm) = e3t_0(:,:,jk) * ( 1._wp + r3t(:,:,Kmm) * tmask(:,:,jk) ) ! Kmm time level e3u(:,:,jk,Kmm) = e3u_0(:,:,jk) * ( 1._wp + r3u(:,:,Kmm) * umask(:,:,jk) ) e3v(:,:,jk,Kmm) = e3v_0(:,:,jk) * ( 1._wp + r3v(:,:,Kmm) * vmask(:,:,jk) ) e3f(:,:,jk) = e3f_0(:,:,jk) * ( 1._wp + r3f(:,:) * fmask(:,:,jk) ) END DO ! DO jk = 1, jpk ! Vertical interpolation of e3t,u,v ! ! The ratio does not have to be masked at w-level e3w (:,:,jk,Kbb) = e3w_0(:,:,jk) * ( 1._wp + r3t(:,:,Kbb) ) ! Kbb time level e3uw(:,:,jk,Kbb) = e3uw_0(:,:,jk) * ( 1._wp + r3u(:,:,Kbb) ) e3vw(:,:,jk,Kbb) = e3vw_0(:,:,jk) * ( 1._wp + r3v(:,:,Kbb) ) e3w (:,:,jk,Kmm) = e3w_0(:,:,jk) * ( 1._wp + r3t(:,:,Kmm) ) ! Kmm time level e3uw(:,:,jk,Kmm) = e3uw_0(:,:,jk) * ( 1._wp + r3u(:,:,Kmm) ) e3vw(:,:,jk,Kmm) = e3vw_0(:,:,jk) * ( 1._wp + r3v(:,:,Kmm) ) END DO ! ! We need to define e3[tuv]_a for AGRIF initialisation (should not be a problem for the restartability...) e3t(:,:,:,Kaa) = e3t(:,:,:,Kmm) e3u(:,:,:,Kaa) = e3u(:,:,:,Kmm) e3v(:,:,:,Kaa) = e3v(:,:,:,Kmm) ! ! !== depth of t and w-point ==! (set the isf depth as it is in the initial timestep) IF( ln_isf ) THEN !** IceShelF cavities ! ! to be created depending of the new names in isf ! ! it should be something like that : (with h_isf = thickness of iceshelf) ! ! in fact currently, h_isf(:,:) is called : risfdep(:,:) !!gm - depth idea 0 : just realize that mask is not needed ===>>>> with ISF, rescale all grid point position below ISF : no mask ! gdept(:,:,1,Kmm) = gdept_0(:,:,1) * ( 1._wp + r3t(:,:,Kmm) ) gdepw(:,:,1,Kmm) = 0._wp ! Initialized to zero one for all gde3w(:,:,1) = gdept(:,:,1,Kmm) - ssh(:,:,Kmm) ! reference to a common level z=0 for hpg DO jk = 2, jpk gdept(:,:,jk,Kmm) = MIN( risfdep(:,:) , gdept_0(:,:,jk) ) & + MAX( 0._wp , gdept_0(:,:,jk)-risfdep(:,:) ) * ( 1._wp + r3t(:,:,Kmm) ) gdepw(:,:,jk,Kmm) = MIN( risfdep(:,:) , gdepw_0(:,:,jk) ) & + MAX( 0._wp , gdepw_0(:,:,jk)-risfdep(:,:) ) * ( 1._wp + r3t(:,:,Kmm) ) gde3w(:,:,jk) = gdept(:,:,jk,Kmm) - ssh(:,:,Kmm) gdept(:,:,jk,Kbb) = MIN( risfdep(:,:) , gdept_0(:,:,jk) ) & + MAX( 0._wp , gdept_0(:,:,jk)-risfdep(:,:) ) * ( 1._wp + r3t(:,:,Kbb) ) gdepw(:,:,jk,Kbb) = MIN( risfdep(:,:) , gdepw_0(:,:,jk) ) & + MAX( 0._wp , gdepw_0(:,:,jk)-risfdep(:,:) ) * ( 1._wp + r3t(:,:,Kbb) ) END DO ! ELSE !** No cavities (all depth rescaled, even inside topography: no mask) ! !!gm idea 0 : just realize that mask is not needed ===>>>> without ISF, rescale all grid point position : no mask ! DO jk = 1, jpk gdept(:,:,jk,Kmm) = gdept_0(:,:,jk) * ( 1._wp + r3t(:,:,Kmm) ) gdepw(:,:,jk,Kmm) = gdepw_0(:,:,jk) * ( 1._wp + r3t(:,:,Kmm) ) gde3w(:,:,jk) = gdept (:,:,jk,Kmm) - ssh(:,:,Kmm) gdept(:,:,jk,Kbb) = gdept_0(:,:,jk) * ( 1._wp + r3t(:,:,Kbb) ) gdepw(:,:,jk,Kbb) = gdepw_0(:,:,jk) * ( 1._wp + r3t(:,:,Kbb) ) END DO ! ENDIF ! ! !== thickness of the water column !! (ocean portion only) ht(:,:) = ht_0(:,:) + ssh(:,:,Kmm) hu(:,:,Kbb) = hu_0(:,:) * ( 1._wp + r3u(:,:,Kbb) ) hu(:,:,Kmm) = hu_0(:,:) * ( 1._wp + r3u(:,:,Kmm) ) hv(:,:,Kbb) = hv_0(:,:) * ( 1._wp + r3v(:,:,Kbb) ) hv(:,:,Kmm) = hv_0(:,:) * ( 1._wp + r3v(:,:,Kmm) ) ! ! !== inverse of water column thickness ==! (u- and v- points) r1_hu(:,:,Kbb) = ssumask(:,:) / ( hu(:,:,Kbb) + 1._wp - ssumask(:,:) ) ! _i mask due to ISF r1_hu(:,:,Kmm) = ssumask(:,:) / ( hu(:,:,Kmm) + 1._wp - ssumask(:,:) ) r1_hv(:,:,Kbb) = ssvmask(:,:) / ( hv(:,:,Kbb) + 1._wp - ssvmask(:,:) ) r1_hv(:,:,Kmm) = ssvmask(:,:) / ( hv(:,:,Kmm) + 1._wp - ssvmask(:,:) ) ! END SUBROUTINE dom_qe_zgr SUBROUTINE dom_qe_sf_nxt( kt, Kbb, Kmm, Kaa, kcall ) !!---------------------------------------------------------------------- !! *** ROUTINE dom_qe_sf_nxt *** !! !! ** Purpose : - compute the after scale factors used in tra_zdf, dynnxt, !! tranxt and dynspg routines !! !! ** Method : - z_star case: Repartition of ssh INCREMENT proportionnaly to the level thickness. !! !! ** Action : - hdiv_lf : restoring towards full baroclinic divergence in z_tilde case !! - tilde_e3t_a: after increment of vertical scale factor !! in z_tilde case !! - e3(t/u/v)_a !! !! Reference : Leclair, M., and Madec, G. 2011, Ocean Modelling. !!---------------------------------------------------------------------- INTEGER, INTENT( in ) :: kt ! time step INTEGER, INTENT( in ) :: Kbb, Kmm, Kaa ! time step INTEGER, INTENT( in ), OPTIONAL :: kcall ! optional argument indicating call sequence ! INTEGER :: ji, jj, jk ! dummy loop indices INTEGER , DIMENSION(3) :: ijk_max, ijk_min ! temporary integers REAL(wp) :: z2dt, z_tmin, z_tmax ! local scalars LOGICAL :: ll_do_bclinic ! local logical REAL(wp), DIMENSION(jpi,jpj) :: zht, z_scale, zwu, zwv, zhdiv !!---------------------------------------------------------------------- ! IF( ln_linssh ) RETURN ! No calculation in linear free surface ! IF( ln_timing ) CALL timing_start('dom_qe_sf_nxt') ! IF( kt == nit000 ) THEN IF(lwp) WRITE(numout,*) IF(lwp) WRITE(numout,*) 'dom_qe_sf_nxt : compute after scale factors' IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~' ENDIF ! ******************************* ! ! After acale factors at t-points ! ! ******************************* ! ! ! --------------------------------------------- ! ! ! z_star coordinate and barotropic z-tilde part ! ! ! --------------------------------------------- ! ! ! ! *********************************** ! ! After scale factors at u- v- points ! ! *********************************** ! ! DO jk = 1, jpkm1 e3t(:,:,jk,Kaa) = e3t_0(:,:,jk) * ( 1._wp + r3t(:,:,Kaa) * tmask(:,:,jk) ) e3u(:,:,jk,Kaa) = e3u_0(:,:,jk) * ( 1._wp + r3u(:,:,Kaa) * umask(:,:,jk) ) e3v(:,:,jk,Kaa) = e3v_0(:,:,jk) * ( 1._wp + r3v(:,:,Kaa) * vmask(:,:,jk) ) END DO ! ! *********************************** ! ! After depths at u- v points ! ! *********************************** ! hu(:,:,Kaa) = hu_0(:,:) * ( 1._wp + r3u(:,:,Kaa) ) hv(:,:,Kaa) = hv_0(:,:) * ( 1._wp + r3v(:,:,Kaa) ) ! ! Inverse of the local depth r1_hu(:,:,Kaa) = ssumask(:,:) / ( hu(:,:,Kaa) + 1._wp - ssumask(:,:) ) r1_hv(:,:,Kaa) = ssvmask(:,:) / ( hv(:,:,Kaa) + 1._wp - ssvmask(:,:) ) ! IF( ln_timing ) CALL timing_stop('dom_qe_sf_nxt') ! END SUBROUTINE dom_qe_sf_nxt SUBROUTINE dom_qe_sf_update( kt, Kbb, Kmm, Kaa ) !!---------------------------------------------------------------------- !! *** ROUTINE dom_qe_sf_update *** !! !! ** Purpose : for z tilde case: compute time filter and swap of scale factors !! compute all depths and related variables for next time step !! write outputs and restart file !! !! ** Method : - reconstruct scale factor at other grid points (interpolate) !! - recompute depths and water height fields !! !! ** Action : - Recompute: !! e3(u/v)_b !! e3w(:,:,:,Kmm) !! e3(u/v)w_b !! e3(u/v)w_n !! gdept(:,:,:,Kmm), gdepw(:,:,:,Kmm) and gde3w !! 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. !!---------------------------------------------------------------------- INTEGER, INTENT( in ) :: kt ! time step INTEGER, INTENT( in ) :: Kbb, Kmm, Kaa ! time level indices ! INTEGER :: ji, jj, jk ! dummy loop indices REAL(wp) :: zcoef ! local scalar !!---------------------------------------------------------------------- ! IF( ln_linssh ) RETURN ! No calculation in linear free surface ! IF( ln_timing ) CALL timing_start('dom_qe_sf_update') ! IF( kt == nit000 ) THEN IF(lwp) WRITE(numout,*) IF(lwp) WRITE(numout,*) 'dom_qe_sf_update : - interpolate scale factors and compute depths for next time step' IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~~~~' ENDIF ! ! Compute all missing vertical scale factor and depths ! ==================================================== ! Horizontal scale factor interpolations ! -------------------------------------- ! - ML - e3u(:,:,:,Kbb) and e3v(:,:,:,Kbb) are already computed in dynnxt ! - JC - hu(:,:,:,Kbb), hv(:,:,:,:,Kbb), hur_b, hvr_b also ! Scale factor computation DO jk = 1, jpk ! Horizontal interpolation e3f(:,:,jk) = e3f_0(:,:,jk) * ( 1._wp + r3f(:,:) * fmask(:,:,jk) ) ! Kmm time level ! ! Vertical interpolation ! ! The ratio does not have to be masked at w-level e3w (:,:,jk,Kmm) = e3w_0(:,:,jk) * ( 1._wp + r3t(:,:,Kmm) ) ! Kmm time level e3uw(:,:,jk,Kmm) = e3uw_0(:,:,jk) * ( 1._wp + r3u(:,:,Kmm) ) e3vw(:,:,jk,Kmm) = e3vw_0(:,:,jk) * ( 1._wp + r3v(:,:,Kmm) ) e3w (:,:,jk,Kbb) = e3w_0(:,:,jk) * ( 1._wp + r3t(:,:,Kbb) ) ! Kbb time level e3uw(:,:,jk,Kbb) = e3uw_0(:,:,jk) * ( 1._wp + r3u(:,:,Kbb) ) e3vw(:,:,jk,Kbb) = e3vw_0(:,:,jk) * ( 1._wp + r3v(:,:,Kbb) ) END DO IF( ln_isf ) THEN !** IceShelF cavities ! ! to be created depending of the new names in isf ! ! it should be something like that : (with h_isf = thickness of iceshelf) ! ! in fact currently, h_isf(:,:) is called : risfdep(:,:) !!gm - depth idea 0 : just realize that mask is not needed ===>>>> with ISF, rescale all grid point position below ISF : no mask ! gdept(:,:,1,Kmm) = gdept_0(:,:,1) * ( 1._wp + r3t(:,:,Kmm) ) gdepw(:,:,1,Kmm) = 0._wp ! Initialized to zero one for all gde3w(:,:,1) = gdept(:,:,1,Kmm) - ssh(:,:,Kmm) ! reference to a common level z=0 for hpg DO jk = 2, jpk gdept(:,:,jk,Kmm) = MIN( risfdep(:,:) , gdept_0(:,:,jk) ) & + MAX( 0._wp , gdept_0(:,:,jk)-risfdep(:,:) ) * ( 1._wp + r3t(:,:,Kmm) ) gdepw(:,:,jk,Kmm) = MIN( risfdep(:,:) , gdepw_0(:,:,jk) ) & + MAX( 0._wp , gdepw_0(:,:,jk)-risfdep(:,:) ) * ( 1._wp + r3t(:,:,Kmm) ) gde3w(:,:,jk) = gdept(:,:,jk,Kmm) - ssh(:,:,Kmm) gdept(:,:,jk,Kbb) = MIN( risfdep(:,:) , gdept_0(:,:,jk) ) & + MAX( 0._wp , gdept_0(:,:,jk)-risfdep(:,:) ) * ( 1._wp + r3t(:,:,Kbb) ) gdepw(:,:,jk,Kbb) = MIN( risfdep(:,:) , gdepw_0(:,:,jk) ) & + MAX( 0._wp , gdepw_0(:,:,jk)-risfdep(:,:) ) * ( 1._wp + r3t(:,:,Kbb) ) END DO ! ELSE !** No cavities (all depth rescaled, even inside topography: no mask) ! !!gm idea 0 : just realize that mask is not needed ===>>>> without ISF, rescale all grid point position : no mask ! DO jk = 1, jpk gdept(:,:,jk,Kmm) = gdept_0(:,:,jk) * ( 1._wp + r3t(:,:,Kmm) ) gdepw(:,:,jk,Kmm) = gdepw_0(:,:,jk) * ( 1._wp + r3t(:,:,Kmm) ) gde3w(:,:,jk) = gdept (:,:,jk,Kmm) - ssh(:,:,Kmm) gdept(:,:,jk,Kbb) = gdept_0(:,:,jk) * ( 1._wp + r3t(:,:,Kbb) ) gdepw(:,:,jk,Kbb) = gdepw_0(:,:,jk) * ( 1._wp + r3t(:,:,Kbb) ) END DO ! ENDIF ! Local depth and Inverse of the local depth of the water ! ------------------------------------------------------- ! ht(:,:) = ht_0(:,:) + ssh(:,:,Kmm) ! write restart file ! ================== IF( lrst_oce ) CALL dom_qe_rst( kt, Kbb, Kmm, 'WRITE' ) ! IF( ln_timing ) CALL timing_stop('dom_qe_sf_update') ! END SUBROUTINE dom_qe_sf_update SUBROUTINE dom_qe_interpol( pe3_in, pe3_out, pout ) !!--------------------------------------------------------------------- !! *** ROUTINE dom_qe_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 !!---------------------------------------------------------------------- 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=*) , INTENT(in ) :: pout ! grid point of out scale factors ! ! = 'U', 'V', 'W, 'F', 'UW' or 'VW' ! INTEGER :: ji, jj, jk ! dummy loop indices REAL(wp) :: zlnwd ! =1./0. when ln_wd_il = T/F !!---------------------------------------------------------------------- ! IF(ln_wd_il) THEN zlnwd = 1.0_wp ELSE zlnwd = 0.0_wp END IF ! SELECT CASE ( pout ) !== type of interpolation ==! ! CASE( 'U' ) !* from T- to U-point : hor. surface weighted mean DO_3D_10_10( 1, jpk ) pe3_out(ji,jj,jk) = 0.5_wp * ( umask(ji,jj,jk) * (1.0_wp - zlnwd) + zlnwd ) * r1_e1e2u(ji,jj) & & * ( e1e2t(ji ,jj) * ( pe3_in(ji ,jj,jk) - e3t_0(ji ,jj,jk) ) & & + e1e2t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) ) END_3D CALL lbc_lnk( 'domqe', pe3_out(:,:,:), 'U', 1._wp ) pe3_out(:,:,:) = pe3_out(:,:,:) + e3u_0(:,:,:) ! CASE( 'V' ) !* from T- to V-point : hor. surface weighted mean DO_3D_10_10( 1, jpk ) pe3_out(ji,jj,jk) = 0.5_wp * ( vmask(ji,jj,jk) * (1.0_wp - zlnwd) + zlnwd ) * r1_e1e2v(ji,jj) & & * ( e1e2t(ji,jj ) * ( pe3_in(ji,jj ,jk) - e3t_0(ji,jj ,jk) ) & & + e1e2t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) ) END_3D CALL lbc_lnk( 'domqe', pe3_out(:,:,:), 'V', 1._wp ) pe3_out(:,:,:) = pe3_out(:,:,:) + e3v_0(:,:,:) ! CASE( 'F' ) !* from U-point to F-point : hor. surface weighted mean DO_3D_10_10( 1, jpk ) pe3_out(ji,jj,jk) = 0.5_wp * ( umask(ji,jj,jk) * umask(ji,jj+1,jk) * (1.0_wp - zlnwd) + zlnwd ) & & * r1_e1e2f(ji,jj) & & * ( e1e2u(ji,jj ) * ( pe3_in(ji,jj ,jk) - e3u_0(ji,jj ,jk) ) & & + e1e2u(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3u_0(ji,jj+1,jk) ) ) END_3D CALL lbc_lnk( 'domqe', pe3_out(:,:,:), 'F', 1._wp ) pe3_out(:,:,:) = pe3_out(:,:,:) + e3f_0(:,:,:) ! CASE( 'W' ) !* from T- to W-point : vertical simple mean ! !zlnwd = 1.0_wp pe3_out(:,:,1) = e3w_0(:,:,1) + pe3_in(:,:,1) - e3t_0(:,:,1) ! - ML - The use of mask in this formulea enables the special treatment of the last w-point without indirect adressing !!gm BUG? use here wmask in case of ISF ? to be checked DO jk = 2, jpk pe3_out(:,:,jk) = e3w_0(:,:,jk) + ( 1.0_wp - 0.5_wp * ( tmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) ) & & * ( pe3_in(:,:,jk-1) - e3t_0(:,:,jk-1) ) & & + 0.5_wp * ( tmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) & & * ( pe3_in(:,:,jk ) - e3t_0(:,:,jk ) ) END DO ! CASE( 'UW' ) !* from U- to UW-point : vertical simple mean ! !zlnwd = 1.0_wp pe3_out(:,:,1) = e3uw_0(:,:,1) + pe3_in(:,:,1) - e3u_0(:,:,1) ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without indirect adressing !!gm BUG? use here wumask in case of ISF ? to be checked DO jk = 2, jpk pe3_out(:,:,jk) = e3uw_0(:,:,jk) + ( 1.0_wp - 0.5_wp * ( umask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) ) & & * ( pe3_in(:,:,jk-1) - e3u_0(:,:,jk-1) ) & & + 0.5_wp * ( umask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) & & * ( pe3_in(:,:,jk ) - e3u_0(:,:,jk ) ) END DO ! CASE( 'VW' ) !* from V- to VW-point : vertical simple mean ! !zlnwd = 1.0_wp pe3_out(:,:,1) = e3vw_0(:,:,1) + pe3_in(:,:,1) - e3v_0(:,:,1) ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without indirect adressing !!gm BUG? use here wvmask in case of ISF ? to be checked DO jk = 2, jpk pe3_out(:,:,jk) = e3vw_0(:,:,jk) + ( 1.0_wp - 0.5_wp * ( vmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) ) & & * ( pe3_in(:,:,jk-1) - e3v_0(:,:,jk-1) ) & & + 0.5_wp * ( vmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) & & * ( pe3_in(:,:,jk ) - e3v_0(:,:,jk ) ) END DO END SELECT ! END SUBROUTINE dom_qe_interpol SUBROUTINE dom_qe_r3c( pssh, pr3t, pr3u, pr3v, pr3f ) !!--------------------------------------------------------------------- !! *** ROUTINE r3c *** !! !! ** Purpose : compute the filtered ratio ssh/h_0 at t-,u-,v-,f-points !! !! ** Method : - compute the ssh at u- and v-points (f-point optional) !! Vector Form : surface weighted averaging !! Flux Form : simple averaging !! - compute the ratio ssh/h_0 at t-,u-,v-pts, (f-pt optional) !!---------------------------------------------------------------------- REAL(wp), DIMENSION(:,:) , INTENT(in ) :: pssh ! sea surface height [m] REAL(wp), DIMENSION(:,:) , INTENT( out) :: pr3t, pr3u, pr3v ! ssh/h0 ratio at t-, u-, v-,points [-] REAL(wp), DIMENSION(:,:), OPTIONAL, INTENT( out) :: pr3f ! ssh/h0 ratio at f-point [-] ! INTEGER :: ji, jj ! dummy loop indices !!---------------------------------------------------------------------- ! ! pr3t(:,:) = pssh(:,:) * r1_ht_0(:,:) !== ratio at t-point ==! ! ! ! !== ratio at u-,v-point ==! ! IF( ln_dynadv_vec ) THEN !- Vector Form (thickness weighted averaging) DO_2D_11_11 pr3u(ji,jj) = 0.5_wp * ( e1e2t(ji ,jj) * pssh(ji ,jj) & & + e1e2t(ji+1,jj) * pssh(ji+1,jj) ) * r1_hu_0(ji,jj) * r1_e1e2u(ji,jj) pr3v(ji,jj) = 0.5_wp * ( e1e2t(ji,jj ) * pssh(ji,jj ) & & + e1e2t(ji,jj+1) * pssh(ji,jj+1) ) * r1_hv_0(ji,jj) * r1_e1e2v(ji,jj) END_2D ELSE !- Flux Form (simple averaging) DO_2D_11_11 pr3u(ji,jj) = 0.5_wp * ( pssh(ji ,jj) + pssh(ji+1,jj) ) * r1_hu_0(ji,jj) pr3v(ji,jj) = 0.5_wp * ( pssh(ji,jj ) + pssh(ji,jj+1) ) * r1_hv_0(ji,jj) END_2D ENDIF ! IF( .NOT.PRESENT( pr3f ) ) THEN !- lbc on ratio at u-, v-points only CALL lbc_lnk_multi( 'dom_qe_r3c', pr3u, 'U', 1._wp, pr3v, 'V', 1._wp ) ! ! ELSE !== ratio at f-point ==! ! IF( ln_dynadv_vec ) THEN !- Vector Form (thickness weighted averaging) DO_2D_01_01 ! start from 1 since lbc_lnk('F') doesn't update the 1st row/line pr3f(ji,jj) = 0.25_wp * ( e1e2t(ji ,jj ) * pssh(ji ,jj ) & & + e1e2t(ji+1,jj ) * pssh(ji+1,jj ) & & + e1e2t(ji ,jj+1) * pssh(ji ,jj+1) & & + e1e2t(ji+1,jj+1) * pssh(ji+1,jj+1) ) * r1_hf_0(ji,jj) * r1_e1e2f(ji,jj) END_2D ELSE !- Flux Form (simple averaging) DO_2D_01_01 ! start from 1 since lbc_lnk('F') doesn't update the 1st row/line pr3f(ji,jj) = 0.25_wp * ( pssh(ji ,jj ) + pssh(ji+1,jj ) & & + pssh(ji ,jj+1) + pssh(ji+1,jj+1) ) * r1_hf_0(ji,jj) END_2D ENDIF ! ! lbc on ratio at u-,v-,f-points CALL lbc_lnk_multi( 'dom_qe_r3c', pr3u, 'U', 1._wp, pr3v, 'V', 1._wp, pr3f, 'F', 1._wp ) ! ENDIF ! END SUBROUTINE dom_qe_r3c SUBROUTINE dom_qe_rst( kt, Kbb, Kmm, cdrw ) !!--------------------------------------------------------------------- !! *** ROUTINE dom_qe_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. !!---------------------------------------------------------------------- INTEGER , INTENT(in) :: kt ! ocean time-step INTEGER , INTENT(in) :: Kbb, Kmm ! ocean time level indices CHARACTER(len=*), INTENT(in) :: cdrw ! "READ"/"WRITE" flag ! INTEGER :: ji, jj, jk INTEGER :: id1, id2 ! local integers !!---------------------------------------------------------------------- ! IF( TRIM(cdrw) == 'READ' ) THEN ! Read/initialise ! ! =============== IF( ln_rstart ) THEN !* Read the restart file CALL rst_read_open ! open the restart file if necessary CALL iom_get( numror, jpdom_autoglo, 'sshn' , ssh(:,:,Kmm), ldxios = lrxios ) ! id1 = iom_varid( numror, 'e3t_b', ldstop = .FALSE. ) id2 = iom_varid( numror, 'e3t_n', ldstop = .FALSE. ) ! ! ! --------- ! ! ! all cases ! ! ! --------- ! ! IF( MIN( id1, id2 ) > 0 ) THEN ! all required arrays exist CALL iom_get( numror, jpdom_autoglo, 'e3t_b', e3t(:,:,:,Kbb), ldxios = lrxios ) CALL iom_get( numror, jpdom_autoglo, 'e3t_n', e3t(:,:,:,Kmm), ldxios = lrxios ) ! needed to restart if land processor not computed IF(lwp) write(numout,*) 'dom_qe_rst : e3t(:,:,:,Kbb) and e3t(:,:,:,Kmm) found in restart files' WHERE ( tmask(:,:,:) == 0.0_wp ) e3t(:,:,:,Kmm) = e3t_0(:,:,:) e3t(:,:,:,Kbb) = e3t_0(:,:,:) END WHERE IF( neuler == 0 ) THEN e3t(:,:,:,Kbb) = e3t(:,:,:,Kmm) ENDIF ELSE IF( id1 > 0 ) THEN IF(lwp) write(numout,*) 'dom_qe_rst WARNING : e3t(:,:,:,Kmm) not found in restart files' IF(lwp) write(numout,*) 'e3t_n set equal to e3t_b.' IF(lwp) write(numout,*) 'neuler is forced to 0' CALL iom_get( numror, jpdom_autoglo, 'e3t_b', e3t(:,:,:,Kbb), ldxios = lrxios ) e3t(:,:,:,Kmm) = e3t(:,:,:,Kbb) neuler = 0 ELSE IF( id2 > 0 ) THEN IF(lwp) write(numout,*) 'dom_qe_rst WARNING : e3t(:,:,:,Kbb) not found in restart files' IF(lwp) write(numout,*) 'e3t_b set equal to e3t_n.' IF(lwp) write(numout,*) 'neuler is forced to 0' CALL iom_get( numror, jpdom_autoglo, 'e3t_n', e3t(:,:,:,Kmm), ldxios = lrxios ) e3t(:,:,:,Kbb) = e3t(:,:,:,Kmm) neuler = 0 ELSE IF(lwp) write(numout,*) 'dom_qe_rst WARNING : e3t(:,:,:,Kmm) not found in restart file' IF(lwp) write(numout,*) 'Compute scale factor from sshn' IF(lwp) write(numout,*) 'neuler is forced to 0' DO jk = 1, jpk e3t(:,:,jk,Kmm) = e3t_0(:,:,jk) * ( ht_0(:,:) + ssh(:,:,Kmm) ) & & / ( ht_0(:,:) + 1._wp - ssmask(:,:) ) * tmask(:,:,jk) & & + e3t_0(:,:,jk) * (1._wp -tmask(:,:,jk)) END DO e3t(:,:,:,Kbb) = e3t(:,:,:,Kmm) neuler = 0 ENDIF ! ELSE !* Initialize at "rest" ! IF( ll_wd ) THEN ! MJB ll_wd edits start here - these are essential ! IF( cn_cfg == 'wad' ) THEN ! Wetting and drying test case CALL usr_def_istate( gdept(:,:,:,Kbb), tmask, ts(:,:,:,:,Kbb), uu(:,:,:,Kbb), vv(:,:,:,Kbb), ssh(:,:,Kbb) ) ts (:,:,:,:,Kmm) = ts (:,:,:,:,Kbb) ! set now values from to before ones ssh (:,:,Kmm) = ssh(:,:,Kbb) uu (:,:,:,Kmm) = uu (:,:,:,Kbb) vv (:,:,:,Kmm) = vv (:,:,:,Kbb) ELSE ! if not test case ssh(:,:,Kmm) = -ssh_ref ssh(:,:,Kbb) = -ssh_ref DO_2D_11_11 IF( ht_0(ji,jj)-ssh_ref < rn_wdmin1 ) THEN ! if total depth is less than min depth ssh(ji,jj,Kbb) = rn_wdmin1 - (ht_0(ji,jj) ) ssh(ji,jj,Kmm) = rn_wdmin1 - (ht_0(ji,jj) ) ENDIF END_2D ENDIF !If test case else ! Adjust vertical metrics for all wad DO jk = 1, jpk e3t(:,:,jk,Kmm) = e3t_0(:,:,jk) * ( 1._wp + r3t(:,:,Kmm) * tmask(:,:,jk) ) END DO e3t(:,:,:,Kbb) = e3t(:,:,:,Kmm) DO ji = 1, jpi DO jj = 1, jpj IF ( ht_0(ji,jj) .LE. 0.0 .AND. NINT( ssmask(ji,jj) ) .EQ. 1) THEN CALL ctl_stop( 'dom_qe_rst: ht_0 must be positive at potentially wet points' ) ENDIF END DO END DO ! ELSE ! ! Just to read set ssh in fact, called latter once vertical grid ! is set up: ! CALL usr_def_istate( gdept_0, tmask, ts(:,:,:,:,Kbb), uu(:,:,:,Kbb), vv(:,:,:,Kbb), ssh(:,:,Kbb) ) ! ! ! DO jk=1,jpk ! e3t(:,:,jk,Kbb) = e3t_0(:,:,jk) * ( ht_0(:,:) + ssh(:,:,Kbb) ) & ! & / ( ht_0(:,:) + 1._wp -ssmask(:,:) ) * tmask(:,:,jk) ! END DO ! e3t(:,:,:,Kmm) = e3t(:,:,:,Kbb) ssh(:,:,Kmm)=0._wp e3t(:,:,:,Kmm)=e3t_0(:,:,:) e3t(:,:,:,Kbb)=e3t_0(:,:,:) ! ENDIF ! end of ll_wd edits ! ENDIF ! ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN ! Create restart file ! ! =================== IF(lwp) WRITE(numout,*) '---- dom_qe_rst ----' IF( lwxios ) CALL iom_swap( cwxios_context ) ! ! --------- ! ! ! all cases ! ! ! --------- ! CALL iom_rstput( kt, nitrst, numrow, 'e3t_b', e3t(:,:,:,Kbb), ldxios = lwxios ) CALL iom_rstput( kt, nitrst, numrow, 'e3t_n', e3t(:,:,:,Kmm), ldxios = lwxios ) ! IF( lwxios ) CALL iom_swap( cxios_context ) ENDIF ! END SUBROUTINE dom_qe_rst SUBROUTINE dom_qe_ctl !!--------------------------------------------------------------------- !! *** ROUTINE dom_qe_ctl *** !! !! ** Purpose : Control the consistency between namelist options !! for vertical coordinate !!---------------------------------------------------------------------- INTEGER :: ioptio, ios !! NAMELIST/nam_vvl/ ln_vvl_zstar, ln_vvl_ztilde, ln_vvl_layer, ln_vvl_ztilde_as_zstar, & & ln_vvl_zstar_at_eqtor , rn_ahe3 , rn_rst_e3t , & & rn_lf_cutoff , rn_zdef_max , ln_vvl_dbg ! not yet implemented: ln_vvl_kepe !!---------------------------------------------------------------------- ! READ ( numnam_ref, nam_vvl, IOSTAT = ios, ERR = 901) 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'nam_vvl in reference namelist' ) READ ( numnam_cfg, nam_vvl, IOSTAT = ios, ERR = 902 ) 902 IF( ios > 0 ) CALL ctl_nam ( ios , 'nam_vvl in configuration namelist' ) IF(lwm) WRITE ( numond, nam_vvl ) ! IF(lwp) THEN ! Namelist print WRITE(numout,*) WRITE(numout,*) 'dom_qe_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,*) ' ztilde as zstar ln_vvl_ztilde_as_zstar = ', ln_vvl_ztilde_as_zstar WRITE(numout,*) ' ztilde near the equator ln_vvl_zstar_at_eqtor = ', ln_vvl_zstar_at_eqtor WRITE(numout,*) ' !' WRITE(numout,*) ' thickness diffusion coefficient rn_ahe3 = ', rn_ahe3 WRITE(numout,*) ' maximum e3t deformation fractional change rn_zdef_max = ', rn_zdef_max IF( ln_vvl_ztilde_as_zstar ) THEN WRITE(numout,*) ' ztilde running in zstar emulation mode (ln_vvl_ztilde_as_zstar=T) ' WRITE(numout,*) ' ignoring namelist timescale parameters and using:' WRITE(numout,*) ' hard-wired : z-tilde to zstar restoration timescale (days)' WRITE(numout,*) ' rn_rst_e3t = 0.e0' WRITE(numout,*) ' hard-wired : z-tilde cutoff frequency of low-pass filter (days)' WRITE(numout,*) ' rn_lf_cutoff = 1.0/rdt' ELSE WRITE(numout,*) ' z-tilde to zstar restoration timescale (days) rn_rst_e3t = ', rn_rst_e3t WRITE(numout,*) ' z-tilde cutoff frequency of low-pass filter (days) rn_lf_cutoff = ', rn_lf_cutoff ENDIF WRITE(numout,*) ' debug prints flag ln_vvl_dbg = ', ln_vvl_dbg ENDIF ! ioptio = 0 ! Parameter control IF( ln_vvl_ztilde_as_zstar ) ln_vvl_ztilde = .true. 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(lwp) THEN ! Print the choice WRITE(numout,*) IF( ln_vvl_zstar ) WRITE(numout,*) ' ==>>> zstar vertical coordinate is used' IF( ln_vvl_ztilde ) WRITE(numout,*) ' ==>>> ztilde vertical coordinate is used' IF( ln_vvl_layer ) WRITE(numout,*) ' ==>>> layer vertical coordinate is used' IF( ln_vvl_ztilde_as_zstar ) WRITE(numout,*) ' ==>>> to emulate a zstar coordinate' ENDIF ! #if defined key_agrif IF( (.NOT.Agrif_Root()).AND.(.NOT.ln_vvl_zstar) ) CALL ctl_stop( 'AGRIF is implemented with zstar coordinate only' ) #endif ! END SUBROUTINE dom_qe_ctl !!====================================================================== END MODULE domqe