New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
domvvl.F90 in NEMO/branches/2019/ENHANCE-02_ISF_nemo/src/OCE/DOM – NEMO

source: NEMO/branches/2019/ENHANCE-02_ISF_nemo/src/OCE/DOM/domvvl.F90 @ 11553

Last change on this file since 11553 was 11553, checked in by mathiot, 5 years ago

ENHANCE-02_ISF: fix coupling issue (ticket #2142)

  • Property svn:keywords set to Id
File size: 55.6 KB
RevLine 
[592]1MODULE domvvl
2   !!======================================================================
3   !!                       ***  MODULE domvvl   ***
4   !! Ocean :
5   !!======================================================================
[1438]6   !! History :  2.0  !  2006-06  (B. Levier, L. Marie)  original code
7   !!            3.1  !  2009-02  (G. Madec, M. Leclair, R. Benshila)  pure z* coordinate
[9019]8   !!            3.3  !  2011-10  (M. Leclair) totally rewrote domvvl: vvl option includes z_star and z_tilde coordinates
[5120]9   !!            3.6  !  2014-11  (P. Mathiot) add ice shelf capability
[592]10   !!----------------------------------------------------------------------
[5836]11
[592]12   !!----------------------------------------------------------------------
[4292]13   !!   dom_vvl_init     : define initial vertical scale factors, depths and column thickness
14   !!   dom_vvl_sf_nxt   : Compute next vertical scale factors
15   !!   dom_vvl_sf_swp   : Swap vertical scale factors and update the vertical grid
16   !!   dom_vvl_interpol : Interpolate vertical scale factors from one grid point to another
17   !!   dom_vvl_rst      : read/write restart file
18   !!   dom_vvl_ctl      : Check the vvl options
19   !!----------------------------------------------------------------------
[592]20   USE oce             ! ocean dynamics and tracers
[11423]21   USE isfcpl
[11521]22   USE isf
[6140]23   USE phycst          ! physical constant
[592]24   USE dom_oce         ! ocean space and time domain
[4292]25   USE sbc_oce         ! ocean surface boundary condition
[6152]26   USE wet_dry         ! wetting and drying
[7646]27   USE usrdef_istate   ! user defined initial state (wad only)
[6140]28   USE restart         ! ocean restart
29   !
[592]30   USE in_out_manager  ! I/O manager
[4292]31   USE iom             ! I/O manager library
[592]32   USE lib_mpp         ! distributed memory computing library
33   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
[3294]34   USE timing          ! Timing
[592]35
36   IMPLICIT NONE
37   PRIVATE
38
[4292]39   PUBLIC  dom_vvl_init       ! called by domain.F90
40   PUBLIC  dom_vvl_sf_nxt     ! called by step.F90
41   PUBLIC  dom_vvl_sf_swp     ! called by step.F90
42   PUBLIC  dom_vvl_interpol   ! called by dynnxt.F90
[592]43
[5836]44   !                                                      !!* Namelist nam_vvl
45   LOGICAL , PUBLIC :: ln_vvl_zstar           = .FALSE.    ! zstar  vertical coordinate
46   LOGICAL , PUBLIC :: ln_vvl_ztilde          = .FALSE.    ! ztilde vertical coordinate
47   LOGICAL , PUBLIC :: ln_vvl_layer           = .FALSE.    ! level  vertical coordinate
48   LOGICAL , PUBLIC :: ln_vvl_ztilde_as_zstar = .FALSE.    ! ztilde vertical coordinate
49   LOGICAL , PUBLIC :: ln_vvl_zstar_at_eqtor  = .FALSE.    ! ztilde vertical coordinate
50   LOGICAL , PUBLIC :: ln_vvl_kepe            = .FALSE.    ! kinetic/potential energy transfer
51   !                                                       ! conservation: not used yet
52   REAL(wp)         :: rn_ahe3                             ! thickness diffusion coefficient
53   REAL(wp)         :: rn_rst_e3t                          ! ztilde to zstar restoration timescale [days]
54   REAL(wp)         :: rn_lf_cutoff                        ! cutoff frequency for low-pass filter  [days]
55   REAL(wp)         :: rn_zdef_max                         ! maximum fractional e3t deformation
56   LOGICAL , PUBLIC :: ln_vvl_dbg = .FALSE.                ! debug control prints
[592]57
[5836]58   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: un_td, vn_td                ! thickness diffusion transport
59   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: hdiv_lf                     ! low frequency part of hz divergence
60   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: tilde_e3t_b, tilde_e3t_n    ! baroclinic scale factors
61   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: tilde_e3t_a, dtilde_e3t_a   ! baroclinic scale factors
62   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   :: frq_rst_e3t                 ! retoring period for scale factors
63   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   :: frq_rst_hdv                 ! retoring period for low freq. divergence
[1438]64
[592]65   !! * Substitutions
66#  include "vectopt_loop_substitute.h90"
67   !!----------------------------------------------------------------------
[10068]68   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[888]69   !! $Id$
[10068]70   !! Software governed by the CeCILL license (see ./LICENSE)
[592]71   !!----------------------------------------------------------------------
[4292]72CONTAINS
73
[2715]74   INTEGER FUNCTION dom_vvl_alloc()
75      !!----------------------------------------------------------------------
[4292]76      !!                ***  FUNCTION dom_vvl_alloc  ***
[2715]77      !!----------------------------------------------------------------------
[6140]78      IF( ln_vvl_zstar )   dom_vvl_alloc = 0
[4292]79      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
[4338]80         ALLOCATE( tilde_e3t_b(jpi,jpj,jpk)  , tilde_e3t_n(jpi,jpj,jpk) , tilde_e3t_a(jpi,jpj,jpk) ,   &
81            &      dtilde_e3t_a(jpi,jpj,jpk) , un_td  (jpi,jpj,jpk)     , vn_td  (jpi,jpj,jpk)     ,   &
82            &      STAT = dom_vvl_alloc        )
[10425]83         CALL mpp_sum ( 'domvvl', dom_vvl_alloc )
84         IF( dom_vvl_alloc /= 0 )   CALL ctl_stop( 'STOP', 'dom_vvl_alloc: failed to allocate arrays' )
[6140]85         un_td = 0._wp
86         vn_td = 0._wp
[4292]87      ENDIF
88      IF( ln_vvl_ztilde ) THEN
89         ALLOCATE( frq_rst_e3t(jpi,jpj) , frq_rst_hdv(jpi,jpj) , hdiv_lf(jpi,jpj,jpk) , STAT= dom_vvl_alloc )
[10425]90         CALL mpp_sum ( 'domvvl', dom_vvl_alloc )
91         IF( dom_vvl_alloc /= 0 )   CALL ctl_stop( 'STOP', 'dom_vvl_alloc: failed to allocate arrays' )
[4292]92      ENDIF
[6140]93      !
[2715]94   END FUNCTION dom_vvl_alloc
95
96
[4292]97   SUBROUTINE dom_vvl_init
[592]98      !!----------------------------------------------------------------------
[4292]99      !!                ***  ROUTINE dom_vvl_init  ***
[592]100      !!                   
[4292]101      !! ** Purpose :  Initialization of all scale factors, depths
102      !!               and water column heights
103      !!
104      !! ** Method  :  - use restart file and/or initialize
105      !!               - interpolate scale factors
106      !!
[6140]107      !! ** Action  : - e3t_(n/b) and tilde_e3t_(n/b)
108      !!              - Regrid: e3(u/v)_n
109      !!                        e3(u/v)_b       
110      !!                        e3w_n           
111      !!                        e3(u/v)w_b     
112      !!                        e3(u/v)w_n     
113      !!                        gdept_n, gdepw_n and gde3w_n
[4292]114      !!              - h(t/u/v)_0
115      !!              - frq_rst_e3t and frq_rst_hdv
116      !!
117      !! Reference  : Leclair, M., and G. Madec, 2011, Ocean Modelling.
[592]118      !!----------------------------------------------------------------------
[6140]119      INTEGER ::   ji, jj, jk
[4292]120      INTEGER ::   ii0, ii1, ij0, ij1
[5120]121      REAL(wp)::   zcoef
[592]122      !!----------------------------------------------------------------------
[6140]123      !
[4292]124      IF(lwp) WRITE(numout,*)
125      IF(lwp) WRITE(numout,*) 'dom_vvl_init : Variable volume activated'
126      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~'
[6140]127      !
128      CALL dom_vvl_ctl     ! choose vertical coordinate (z_star, z_tilde or layer)
129      !
130      !                    ! Allocate module arrays
[4292]131      IF( dom_vvl_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dom_vvl_init : unable to allocate arrays' )
[6140]132      !
133      !                    ! Read or initialize e3t_(b/n), tilde_e3t_(b/n) and hdiv_lf
[4292]134      CALL dom_vvl_rst( nit000, 'READ' )
[7753]135      e3t_a(:,:,jpk) = e3t_0(:,:,jpk)  ! last level always inside the sea floor set one for all
[6140]136      !
137      !                    !== Set of all other vertical scale factors  ==!  (now and before)
138      !                                ! Horizontal interpolation of e3t
139      CALL dom_vvl_interpol( e3t_b(:,:,:), e3u_b(:,:,:), 'U' )    ! from T to U
140      CALL dom_vvl_interpol( e3t_n(:,:,:), e3u_n(:,:,:), 'U' )
141      CALL dom_vvl_interpol( e3t_b(:,:,:), e3v_b(:,:,:), 'V' )    ! from T to V
142      CALL dom_vvl_interpol( e3t_n(:,:,:), e3v_n(:,:,:), 'V' )
143      CALL dom_vvl_interpol( e3u_n(:,:,:), e3f_n(:,:,:), 'F' )    ! from U to F
144      !                                ! Vertical interpolation of e3t,u,v
145      CALL dom_vvl_interpol( e3t_n(:,:,:), e3w_n (:,:,:), 'W'  )  ! from T to W
146      CALL dom_vvl_interpol( e3t_b(:,:,:), e3w_b (:,:,:), 'W'  )
147      CALL dom_vvl_interpol( e3u_n(:,:,:), e3uw_n(:,:,:), 'UW' )  ! from U to UW
148      CALL dom_vvl_interpol( e3u_b(:,:,:), e3uw_b(:,:,:), 'UW' )
149      CALL dom_vvl_interpol( e3v_n(:,:,:), e3vw_n(:,:,:), 'VW' )  ! from V to UW
150      CALL dom_vvl_interpol( e3v_b(:,:,:), e3vw_b(:,:,:), 'VW' )
[9449]151
152      ! We need to define e3[tuv]_a for AGRIF initialisation (should not be a problem for the restartability...)
153      e3t_a(:,:,:) = e3t_n(:,:,:)
154      e3u_a(:,:,:) = e3u_n(:,:,:)
155      e3v_a(:,:,:) = e3v_n(:,:,:)
[6140]156      !
157      !                    !==  depth of t and w-point  ==!   (set the isf depth as it is in the initial timestep)
[7753]158      gdept_n(:,:,1) = 0.5_wp * e3w_n(:,:,1)       ! reference to the ocean surface (used for MLD and light penetration)
159      gdepw_n(:,:,1) = 0.0_wp
160      gde3w_n(:,:,1) = gdept_n(:,:,1) - sshn(:,:)  ! reference to a common level z=0 for hpg
161      gdept_b(:,:,1) = 0.5_wp * e3w_b(:,:,1)
162      gdepw_b(:,:,1) = 0.0_wp
[6140]163      DO jk = 2, jpk                               ! vertical sum
[5120]164         DO jj = 1,jpj
165            DO ji = 1,jpi
[6140]166               !    zcoef = tmask - wmask    ! 0 everywhere tmask = wmask, ie everywhere expect at jk = mikt
167               !                             ! 1 everywhere from mbkt to mikt + 1 or 1 (if no isf)
168               !                             ! 0.5 where jk = mikt     
169!!gm ???????   BUG ?  gdept_n as well as gde3w_n  does not include the thickness of ISF ??
170               zcoef = ( tmask(ji,jj,jk) - wmask(ji,jj,jk) )
171               gdepw_n(ji,jj,jk) = gdepw_n(ji,jj,jk-1) + e3t_n(ji,jj,jk-1)
172               gdept_n(ji,jj,jk) =      zcoef  * ( gdepw_n(ji,jj,jk  ) + 0.5 * e3w_n(ji,jj,jk))  &
173                  &                + (1-zcoef) * ( gdept_n(ji,jj,jk-1) +       e3w_n(ji,jj,jk)) 
174               gde3w_n(ji,jj,jk) = gdept_n(ji,jj,jk) - sshn(ji,jj)
175               gdepw_b(ji,jj,jk) = gdepw_b(ji,jj,jk-1) + e3t_b(ji,jj,jk-1)
176               gdept_b(ji,jj,jk) =      zcoef  * ( gdepw_b(ji,jj,jk  ) + 0.5 * e3w_b(ji,jj,jk))  &
177                  &                + (1-zcoef) * ( gdept_b(ji,jj,jk-1) +       e3w_b(ji,jj,jk)) 
[4990]178            END DO
179         END DO
[592]180      END DO
[6140]181      !
182      !                    !==  thickness of the water column  !!   (ocean portion only)
[7753]183      ht_n(:,:) = e3t_n(:,:,1) * tmask(:,:,1)   !!gm  BUG  :  this should be 1/2 * e3w(k=1) ....
184      hu_b(:,:) = e3u_b(:,:,1) * umask(:,:,1)
185      hu_n(:,:) = e3u_n(:,:,1) * umask(:,:,1)
186      hv_b(:,:) = e3v_b(:,:,1) * vmask(:,:,1)
187      hv_n(:,:) = e3v_n(:,:,1) * vmask(:,:,1)
[6140]188      DO jk = 2, jpkm1
[7753]189         ht_n(:,:) = ht_n(:,:) + e3t_n(:,:,jk) * tmask(:,:,jk)
190         hu_b(:,:) = hu_b(:,:) + e3u_b(:,:,jk) * umask(:,:,jk)
191         hu_n(:,:) = hu_n(:,:) + e3u_n(:,:,jk) * umask(:,:,jk)
192         hv_b(:,:) = hv_b(:,:) + e3v_b(:,:,jk) * vmask(:,:,jk)
193         hv_n(:,:) = hv_n(:,:) + e3v_n(:,:,jk) * vmask(:,:,jk)
[4370]194      END DO
[6140]195      !
196      !                    !==  inverse of water column thickness   ==!   (u- and v- points)
[7753]197      r1_hu_b(:,:) = ssumask(:,:) / ( hu_b(:,:) + 1._wp - ssumask(:,:) )    ! _i mask due to ISF
198      r1_hu_n(:,:) = ssumask(:,:) / ( hu_n(:,:) + 1._wp - ssumask(:,:) )
199      r1_hv_b(:,:) = ssvmask(:,:) / ( hv_b(:,:) + 1._wp - ssvmask(:,:) )
200      r1_hv_n(:,:) = ssvmask(:,:) / ( hv_n(:,:) + 1._wp - ssvmask(:,:) )
201
[6140]202      !                    !==   z_tilde coordinate case  ==!   (Restoring frequencies)
[4292]203      IF( ln_vvl_ztilde ) THEN
[6140]204!!gm : idea: add here a READ in a file of custumized restoring frequency
205         !                                   ! Values in days provided via the namelist
206         !                                   ! use rsmall to avoid possible division by zero errors with faulty settings
[7753]207         frq_rst_e3t(:,:) = 2._wp * rpi / ( MAX( rn_rst_e3t  , rsmall ) * 86400.0_wp )
208         frq_rst_hdv(:,:) = 2._wp * rpi / ( MAX( rn_lf_cutoff, rsmall ) * 86400.0_wp )
[6140]209         !
210         IF( ln_vvl_ztilde_as_zstar ) THEN   ! z-star emulation using z-tile
[7753]211            frq_rst_e3t(:,:) = 0._wp               !Ignore namelist settings
212            frq_rst_hdv(:,:) = 1._wp / rdt
[4292]213         ENDIF
[6140]214         IF ( ln_vvl_zstar_at_eqtor ) THEN   ! use z-star in vicinity of the Equator
[4292]215            DO jj = 1, jpj
216               DO ji = 1, jpi
[6140]217!!gm  case |gphi| >= 6 degrees is useless   initialized just above by default
[4292]218                  IF( ABS(gphit(ji,jj)) >= 6.) THEN
219                     ! values outside the equatorial band and transition zone (ztilde)
220                     frq_rst_e3t(ji,jj) =  2.0_wp * rpi / ( MAX( rn_rst_e3t  , rsmall ) * 86400.e0_wp )
221                     frq_rst_hdv(ji,jj) =  2.0_wp * rpi / ( MAX( rn_lf_cutoff, rsmall ) * 86400.e0_wp )
[6140]222                  ELSEIF( ABS(gphit(ji,jj)) <= 2.5) THEN    ! Equator strip ==> z-star
[4292]223                     ! values inside the equatorial band (ztilde as zstar)
224                     frq_rst_e3t(ji,jj) =  0.0_wp
225                     frq_rst_hdv(ji,jj) =  1.0_wp / rdt
[6140]226                  ELSE                                      ! transition band (2.5 to 6 degrees N/S)
227                     !                                      ! (linearly transition from z-tilde to z-star)
[4292]228                     frq_rst_e3t(ji,jj) = 0.0_wp + (frq_rst_e3t(ji,jj)-0.0_wp)*0.5_wp   &
229                        &            * (  1.0_wp - COS( rad*(ABS(gphit(ji,jj))-2.5_wp)  &
230                        &                                          * 180._wp / 3.5_wp ) )
231                     frq_rst_hdv(ji,jj) = (1.0_wp / rdt)                                &
232                        &            + (  frq_rst_hdv(ji,jj)-(1.e0_wp / rdt) )*0.5_wp   &
233                        &            * (  1._wp  - COS( rad*(ABS(gphit(ji,jj))-2.5_wp)  &
234                        &                                          * 180._wp / 3.5_wp ) )
235                  ENDIF
236               END DO
237            END DO
[10213]238            IF( cn_cfg == "orca" .OR. cn_cfg == "ORCA" ) THEN
239               IF( nn_cfg == 3 ) THEN   ! ORCA2: Suppress ztilde in the Foxe Basin for ORCA2
240                  ii0 = 103   ;   ii1 = 111       
241                  ij0 = 128   ;   ij1 = 135   ;   
242                  frq_rst_e3t( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) =  0.0_wp
243                  frq_rst_hdv( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) =  1.e0_wp / rdt
244               ENDIF
[4292]245            ENDIF
246         ENDIF
247      ENDIF
[6140]248      !
[9367]249      IF(lwxios) THEN
250! define variables in restart file when writing with XIOS
251         CALL iom_set_rstw_var_active('e3t_b')
252         CALL iom_set_rstw_var_active('e3t_n')
253         !                                           ! ----------------------- !
254         IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN  ! z_tilde and layer cases !
255            !                                        ! ----------------------- !
256            CALL iom_set_rstw_var_active('tilde_e3t_b')
257            CALL iom_set_rstw_var_active('tilde_e3t_n')
258         END IF
259         !                                           ! -------------!   
260         IF( ln_vvl_ztilde ) THEN                    ! z_tilde case !
261            !                                        ! ------------ !
262            CALL iom_set_rstw_var_active('hdiv_lf')
263         ENDIF
264         !
265      ENDIF
266      !
[4292]267   END SUBROUTINE dom_vvl_init
268
269
[4338]270   SUBROUTINE dom_vvl_sf_nxt( kt, kcall ) 
[4292]271      !!----------------------------------------------------------------------
272      !!                ***  ROUTINE dom_vvl_sf_nxt  ***
273      !!                   
274      !! ** Purpose :  - compute the after scale factors used in tra_zdf, dynnxt,
275      !!                 tranxt and dynspg routines
276      !!
277      !! ** Method  :  - z_star case:  Repartition of ssh INCREMENT proportionnaly to the level thickness.
278      !!               - z_tilde_case: after scale factor increment =
279      !!                                    high frequency part of horizontal divergence
280      !!                                  + retsoring towards the background grid
281      !!                                  + thickness difusion
282      !!                               Then repartition of ssh INCREMENT proportionnaly
283      !!                               to the "baroclinic" level thickness.
284      !!
285      !! ** Action  :  - hdiv_lf    : restoring towards full baroclinic divergence in z_tilde case
286      !!               - tilde_e3t_a: after increment of vertical scale factor
287      !!                              in z_tilde case
[6140]288      !!               - e3(t/u/v)_a
[4292]289      !!
290      !! Reference  : Leclair, M., and Madec, G. 2011, Ocean Modelling.
291      !!----------------------------------------------------------------------
[6140]292      INTEGER, INTENT( in )           ::   kt      ! time step
293      INTEGER, INTENT( in ), OPTIONAL ::   kcall   ! optional argument indicating call sequence
294      !
295      INTEGER                ::   ji, jj, jk            ! dummy loop indices
296      INTEGER , DIMENSION(3) ::   ijk_max, ijk_min      ! temporary integers
297      REAL(wp)               ::   z2dt, z_tmin, z_tmax  ! local scalars
298      LOGICAL                ::   ll_do_bclinic         ! local logical
[9019]299      REAL(wp), DIMENSION(jpi,jpj)     ::   zht, z_scale, zwu, zwv, zhdiv
300      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   ze3t
[4292]301      !!----------------------------------------------------------------------
[6140]302      !
303      IF( ln_linssh )   RETURN      ! No calculation in linear free surface
304      !
[9019]305      IF( ln_timing )   CALL timing_start('dom_vvl_sf_nxt')
[6140]306      !
307      IF( kt == nit000 ) THEN
[4292]308         IF(lwp) WRITE(numout,*)
309         IF(lwp) WRITE(numout,*) 'dom_vvl_sf_nxt : compute after scale factors'
310         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~'
311      ENDIF
312
[4338]313      ll_do_bclinic = .TRUE.
314      IF( PRESENT(kcall) ) THEN
[6140]315         IF( kcall == 2 .AND. ln_vvl_ztilde )   ll_do_bclinic = .FALSE.
[4338]316      ENDIF
317
[4292]318      ! ******************************* !
319      ! After acale factors at t-points !
320      ! ******************************* !
[4338]321      !                                                ! --------------------------------------------- !
[6140]322      !                                                ! z_star coordinate and barotropic z-tilde part !
[4338]323      !                                                ! --------------------------------------------- !
[6140]324      !
[7753]325      z_scale(:,:) = ( ssha(:,:) - sshb(:,:) ) * ssmask(:,:) / ( ht_0(:,:) + sshn(:,:) + 1. - ssmask(:,:) )
[4338]326      DO jk = 1, jpkm1
[7753]327         ! formally this is the same as e3t_a = e3t_0*(1+ssha/ht_0)
328         e3t_a(:,:,jk) = e3t_b(:,:,jk) + e3t_n(:,:,jk) * z_scale(:,:) * tmask(:,:,jk)
[4338]329      END DO
[6140]330      !
[4338]331      IF( ln_vvl_ztilde .OR. ln_vvl_layer .AND. ll_do_bclinic ) THEN   ! z_tilde or layer coordinate !
332         !                                                            ! ------baroclinic part------ !
[4292]333         ! I - initialization
334         ! ==================
335
336         ! 1 - barotropic divergence
337         ! -------------------------
[7753]338         zhdiv(:,:) = 0._wp
339         zht(:,:)   = 0._wp
[4292]340         DO jk = 1, jpkm1
[7753]341            zhdiv(:,:) = zhdiv(:,:) + e3t_n(:,:,jk) * hdivn(:,:,jk)
342            zht  (:,:) = zht  (:,:) + e3t_n(:,:,jk) * tmask(:,:,jk)
[592]343         END DO
[7753]344         zhdiv(:,:) = zhdiv(:,:) / ( zht(:,:) + 1. - tmask_i(:,:) )
[2528]345
[4292]346         ! 2 - Low frequency baroclinic horizontal divergence  (z-tilde case only)
347         ! --------------------------------------------------
348         IF( ln_vvl_ztilde ) THEN
[6140]349            IF( kt > nit000 ) THEN
[4292]350               DO jk = 1, jpkm1
[7753]351                  hdiv_lf(:,:,jk) = hdiv_lf(:,:,jk) - rdt * frq_rst_hdv(:,:)   &
352                     &          * ( hdiv_lf(:,:,jk) - e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) )
[4292]353               END DO
354            ENDIF
[6140]355         ENDIF
[3294]356
[4292]357         ! II - after z_tilde increments of vertical scale factors
358         ! =======================================================
[7753]359         tilde_e3t_a(:,:,:) = 0._wp  ! tilde_e3t_a used to store tendency terms
[4292]360
361         ! 1 - High frequency divergence term
362         ! ----------------------------------
363         IF( ln_vvl_ztilde ) THEN     ! z_tilde case
364            DO jk = 1, jpkm1
[7753]365               tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - ( e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) - hdiv_lf(:,:,jk) )
[4292]366            END DO
367         ELSE                         ! layer case
368            DO jk = 1, jpkm1
[7753]369               tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) -   e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) * tmask(:,:,jk)
[4292]370            END DO
[6140]371         ENDIF
[4292]372
373         ! 2 - Restoring term (z-tilde case only)
374         ! ------------------
375         IF( ln_vvl_ztilde ) THEN
376            DO jk = 1, jpk
[7753]377               tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - frq_rst_e3t(:,:) * tilde_e3t_b(:,:,jk)
[4292]378            END DO
[6140]379         ENDIF
[4292]380
381         ! 3 - Thickness diffusion term
382         ! ----------------------------
[7753]383         zwu(:,:) = 0._wp
384         zwv(:,:) = 0._wp
[6140]385         DO jk = 1, jpkm1        ! a - first derivative: diffusive fluxes
[4292]386            DO jj = 1, jpjm1
387               DO ji = 1, fs_jpim1   ! vector opt.
[5836]388                  un_td(ji,jj,jk) = rn_ahe3 * umask(ji,jj,jk) * e2_e1u(ji,jj)           &
389                     &            * ( tilde_e3t_b(ji,jj,jk) - tilde_e3t_b(ji+1,jj  ,jk) )
390                  vn_td(ji,jj,jk) = rn_ahe3 * vmask(ji,jj,jk) * e1_e2v(ji,jj)           & 
391                     &            * ( tilde_e3t_b(ji,jj,jk) - tilde_e3t_b(ji  ,jj+1,jk) )
[4292]392                  zwu(ji,jj) = zwu(ji,jj) + un_td(ji,jj,jk)
393                  zwv(ji,jj) = zwv(ji,jj) + vn_td(ji,jj,jk)
394               END DO
395            END DO
396         END DO
[6140]397         DO jj = 1, jpj          ! b - correction for last oceanic u-v points
[4292]398            DO ji = 1, jpi
399               un_td(ji,jj,mbku(ji,jj)) = un_td(ji,jj,mbku(ji,jj)) - zwu(ji,jj)
400               vn_td(ji,jj,mbkv(ji,jj)) = vn_td(ji,jj,mbkv(ji,jj)) - zwv(ji,jj)
401            END DO
402         END DO
[6140]403         DO jk = 1, jpkm1        ! c - second derivative: divergence of diffusive fluxes
[4292]404            DO jj = 2, jpjm1
405               DO ji = fs_2, fs_jpim1   ! vector opt.
406                  tilde_e3t_a(ji,jj,jk) = tilde_e3t_a(ji,jj,jk) + (   un_td(ji-1,jj  ,jk) - un_td(ji,jj,jk)    &
407                     &                                          +     vn_td(ji  ,jj-1,jk) - vn_td(ji,jj,jk)    &
[5836]408                     &                                            ) * r1_e1e2t(ji,jj)
[4292]409               END DO
410            END DO
411         END DO
[6140]412         !                       ! d - thickness diffusion transport: boundary conditions
413         !                             (stored for tracer advction and continuity equation)
[10425]414         CALL lbc_lnk_multi( 'domvvl', un_td , 'U' , -1._wp, vn_td , 'V' , -1._wp)
[4292]415
416         ! 4 - Time stepping of baroclinic scale factors
417         ! ---------------------------------------------
418         ! Leapfrog time stepping
419         ! ~~~~~~~~~~~~~~~~~~~~~~
420         IF( neuler == 0 .AND. kt == nit000 ) THEN
421            z2dt =  rdt
422         ELSE
423            z2dt = 2.0_wp * rdt
424         ENDIF
[10425]425         CALL lbc_lnk( 'domvvl', tilde_e3t_a(:,:,:), 'T', 1._wp )
[7753]426         tilde_e3t_a(:,:,:) = tilde_e3t_b(:,:,:) + z2dt * tmask(:,:,:) * tilde_e3t_a(:,:,:)
[4292]427
428         ! Maximum deformation control
429         ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~
[7753]430         ze3t(:,:,jpk) = 0._wp
[4292]431         DO jk = 1, jpkm1
[7753]432            ze3t(:,:,jk) = tilde_e3t_a(:,:,jk) / e3t_0(:,:,jk) * tmask(:,:,jk) * tmask_i(:,:)
[4292]433         END DO
434         z_tmax = MAXVAL( ze3t(:,:,:) )
[10425]435         CALL mpp_max( 'domvvl', z_tmax )                 ! max over the global domain
[4292]436         z_tmin = MINVAL( ze3t(:,:,:) )
[10425]437         CALL mpp_min( 'domvvl', z_tmin )                 ! min over the global domain
[4292]438         ! - ML - test: for the moment, stop simulation for too large e3_t variations
[6140]439         IF( ( z_tmax >  rn_zdef_max ) .OR. ( z_tmin < - rn_zdef_max ) ) THEN
[4292]440            IF( lk_mpp ) THEN
[10425]441               CALL mpp_maxloc( 'domvvl', ze3t, tmask, z_tmax, ijk_max )
442               CALL mpp_minloc( 'domvvl', ze3t, tmask, z_tmin, ijk_min )
[4292]443            ELSE
444               ijk_max = MAXLOC( ze3t(:,:,:) )
445               ijk_max(1) = ijk_max(1) + nimpp - 1
446               ijk_max(2) = ijk_max(2) + njmpp - 1
447               ijk_min = MINLOC( ze3t(:,:,:) )
448               ijk_min(1) = ijk_min(1) + nimpp - 1
449               ijk_min(2) = ijk_min(2) + njmpp - 1
450            ENDIF
451            IF (lwp) THEN
452               WRITE(numout, *) 'MAX( tilde_e3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmax
453               WRITE(numout, *) 'at i, j, k=', ijk_max
454               WRITE(numout, *) 'MIN( tilde_e3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmin
455               WRITE(numout, *) 'at i, j, k=', ijk_min           
[10425]456               CALL ctl_stop( 'STOP', 'MAX( ABS( tilde_e3t_a(:,:,: ) ) / e3t_0(:,:,:) ) too high')
[4292]457            ENDIF
458         ENDIF
459         ! - ML - end test
460         ! - ML - Imposing these limits will cause a baroclinicity error which is corrected for below
[7753]461         tilde_e3t_a(:,:,:) = MIN( tilde_e3t_a(:,:,:),   rn_zdef_max * e3t_0(:,:,:) )
462         tilde_e3t_a(:,:,:) = MAX( tilde_e3t_a(:,:,:), - rn_zdef_max * e3t_0(:,:,:) )
[4292]463
[4338]464         !
465         ! "tilda" change in the after scale factor
[4292]466         ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
[4338]467         DO jk = 1, jpkm1
[7753]468            dtilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - tilde_e3t_b(:,:,jk)
[4338]469         END DO
[4292]470         ! III - Barotropic repartition of the sea surface height over the baroclinic profile
471         ! ==================================================================================
[4338]472         ! add ( ssh increment + "baroclinicity error" ) proportionly to e3t(n)
[4292]473         ! - ML - baroclinicity error should be better treated in the future
474         !        i.e. locally and not spread over the water column.
475         !        (keep in mind that the idea is to reduce Eulerian velocity as much as possible)
[7753]476         zht(:,:) = 0.
[4292]477         DO jk = 1, jpkm1
[7753]478            zht(:,:)  = zht(:,:) + tilde_e3t_a(:,:,jk) * tmask(:,:,jk)
[4292]479         END DO
[7753]480         z_scale(:,:) =  - zht(:,:) / ( ht_0(:,:) + sshn(:,:) + 1. - ssmask(:,:) )
[4292]481         DO jk = 1, jpkm1
[7753]482            dtilde_e3t_a(:,:,jk) = dtilde_e3t_a(:,:,jk) + e3t_n(:,:,jk) * z_scale(:,:) * tmask(:,:,jk)
[4292]483         END DO
[7753]484
[4292]485      ENDIF
486
[4338]487      IF( ln_vvl_ztilde .OR. ln_vvl_layer )  THEN   ! z_tilde or layer coordinate !
488      !                                           ! ---baroclinic part--------- !
489         DO jk = 1, jpkm1
[7753]490            e3t_a(:,:,jk) = e3t_a(:,:,jk) + dtilde_e3t_a(:,:,jk) * tmask(:,:,jk)
[4338]491         END DO
492      ENDIF
493
494      IF( ln_vvl_dbg .AND. .NOT. ll_do_bclinic ) THEN   ! - ML - test: control prints for debuging
[4292]495         !
496         IF( lwp ) WRITE(numout, *) 'kt =', kt
497         IF ( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
498            z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( zht(:,:) ) )
[10425]499            CALL mpp_max( 'domvvl', z_tmax )                             ! max over the global domain
[4292]500            IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(SUM(tilde_e3t_a))) =', z_tmax
501         END IF
502         !
[7753]503         zht(:,:) = 0.0_wp
[4292]504         DO jk = 1, jpkm1
[7753]505            zht(:,:) = zht(:,:) + e3t_n(:,:,jk) * tmask(:,:,jk)
[4292]506         END DO
507         z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + sshn(:,:) - zht(:,:) ) )
[10425]508         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
[6140]509         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+sshn-SUM(e3t_n))) =', z_tmax
[4292]510         !
[7753]511         zht(:,:) = 0.0_wp
[4292]512         DO jk = 1, jpkm1
[7753]513            zht(:,:) = zht(:,:) + e3t_a(:,:,jk) * tmask(:,:,jk)
[4292]514         END DO
515         z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + ssha(:,:) - zht(:,:) ) )
[10425]516         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
[6140]517         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+ssha-SUM(e3t_a))) =', z_tmax
[4292]518         !
[7753]519         zht(:,:) = 0.0_wp
[4292]520         DO jk = 1, jpkm1
[7753]521            zht(:,:) = zht(:,:) + e3t_b(:,:,jk) * tmask(:,:,jk)
[4292]522         END DO
523         z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + sshb(:,:) - zht(:,:) ) )
[10425]524         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
[6140]525         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+sshb-SUM(e3t_b))) =', z_tmax
[4292]526         !
527         z_tmax = MAXVAL( tmask(:,:,1) *  ABS( sshb(:,:) ) )
[10425]528         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
[4292]529         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(sshb))) =', z_tmax
530         !
531         z_tmax = MAXVAL( tmask(:,:,1) *  ABS( sshn(:,:) ) )
[10425]532         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
[4292]533         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(sshn))) =', z_tmax
534         !
535         z_tmax = MAXVAL( tmask(:,:,1) *  ABS( ssha(:,:) ) )
[10425]536         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
[4292]537         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ssha))) =', z_tmax
538      END IF
539
540      ! *********************************** !
541      ! After scale factors at u- v- points !
542      ! *********************************** !
543
[6140]544      CALL dom_vvl_interpol( e3t_a(:,:,:), e3u_a(:,:,:), 'U' )
545      CALL dom_vvl_interpol( e3t_a(:,:,:), e3v_a(:,:,:), 'V' )
[4292]546
[4370]547      ! *********************************** !
548      ! After depths at u- v points         !
549      ! *********************************** !
550
[7753]551      hu_a(:,:) = e3u_a(:,:,1) * umask(:,:,1)
552      hv_a(:,:) = e3v_a(:,:,1) * vmask(:,:,1)
[6140]553      DO jk = 2, jpkm1
[7753]554         hu_a(:,:) = hu_a(:,:) + e3u_a(:,:,jk) * umask(:,:,jk)
555         hv_a(:,:) = hv_a(:,:) + e3v_a(:,:,jk) * vmask(:,:,jk)
[4370]556      END DO
557      !                                        ! Inverse of the local depth
[6140]558!!gm BUG ?  don't understand the use of umask_i here .....
[7753]559      r1_hu_a(:,:) = ssumask(:,:) / ( hu_a(:,:) + 1._wp - ssumask(:,:) )
560      r1_hv_a(:,:) = ssvmask(:,:) / ( hv_a(:,:) + 1._wp - ssvmask(:,:) )
[6140]561      !
[9019]562      IF( ln_timing )   CALL timing_stop('dom_vvl_sf_nxt')
[6140]563      !
[4292]564   END SUBROUTINE dom_vvl_sf_nxt
565
566
567   SUBROUTINE dom_vvl_sf_swp( kt )
[3294]568      !!----------------------------------------------------------------------
[4292]569      !!                ***  ROUTINE dom_vvl_sf_swp  ***
[3294]570      !!                   
[4292]571      !! ** Purpose :  compute time filter and swap of scale factors
572      !!               compute all depths and related variables for next time step
573      !!               write outputs and restart file
[3294]574      !!
[4292]575      !! ** Method  :  - swap of e3t with trick for volume/tracer conservation
576      !!               - reconstruct scale factor at other grid points (interpolate)
577      !!               - recompute depths and water height fields
578      !!
[6140]579      !! ** Action  :  - e3t_(b/n), tilde_e3t_(b/n) and e3(u/v)_n ready for next time step
[4292]580      !!               - Recompute:
[6140]581      !!                    e3(u/v)_b       
582      !!                    e3w_n           
583      !!                    e3(u/v)w_b     
584      !!                    e3(u/v)w_n     
585      !!                    gdept_n, gdepw_n  and gde3w_n
[4292]586      !!                    h(u/v) and h(u/v)r
587      !!
588      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
589      !!              Leclair, M., and G. Madec, 2011, Ocean Modelling.
[3294]590      !!----------------------------------------------------------------------
[6140]591      INTEGER, INTENT( in ) ::   kt   ! time step
592      !
593      INTEGER  ::   ji, jj, jk   ! dummy loop indices
594      REAL(wp) ::   zcoef        ! local scalar
[3294]595      !!----------------------------------------------------------------------
[6140]596      !
597      IF( ln_linssh )   RETURN      ! No calculation in linear free surface
598      !
[9019]599      IF( ln_timing )   CALL timing_start('dom_vvl_sf_swp')
[3294]600      !
[4292]601      IF( kt == nit000 )   THEN
602         IF(lwp) WRITE(numout,*)
603         IF(lwp) WRITE(numout,*) 'dom_vvl_sf_swp : - time filter and swap of scale factors'
604         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~   - interpolate scale factors and compute depths for next time step'
[3294]605      ENDIF
[4292]606      !
607      ! Time filter and swap of scale factors
608      ! =====================================
[6140]609      ! - ML - e3(t/u/v)_b are allready computed in dynnxt.
[4292]610      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
611         IF( neuler == 0 .AND. kt == nit000 ) THEN
[7753]612            tilde_e3t_b(:,:,:) = tilde_e3t_n(:,:,:)
[4292]613         ELSE
[7753]614            tilde_e3t_b(:,:,:) = tilde_e3t_n(:,:,:) & 
615            &         + atfp * ( tilde_e3t_b(:,:,:) - 2.0_wp * tilde_e3t_n(:,:,:) + tilde_e3t_a(:,:,:) )
[4292]616         ENDIF
[7753]617         tilde_e3t_n(:,:,:) = tilde_e3t_a(:,:,:)
[4292]618      ENDIF
[7753]619      gdept_b(:,:,:) = gdept_n(:,:,:)
620      gdepw_b(:,:,:) = gdepw_n(:,:,:)
[4488]621
[7753]622      e3t_n(:,:,:) = e3t_a(:,:,:)
623      e3u_n(:,:,:) = e3u_a(:,:,:)
624      e3v_n(:,:,:) = e3v_a(:,:,:)
625
[4292]626      ! Compute all missing vertical scale factor and depths
627      ! ====================================================
628      ! Horizontal scale factor interpolations
629      ! --------------------------------------
[6140]630      ! - ML - e3u_b and e3v_b are allready computed in dynnxt
[4370]631      ! - JC - hu_b, hv_b, hur_b, hvr_b also
[6140]632     
633      CALL dom_vvl_interpol( e3u_n(:,:,:), e3f_n(:,:,:), 'F'  )
634     
[4292]635      ! Vertical scale factor interpolations
[6140]636      CALL dom_vvl_interpol( e3t_n(:,:,:),  e3w_n(:,:,:), 'W'  )
637      CALL dom_vvl_interpol( e3u_n(:,:,:), e3uw_n(:,:,:), 'UW' )
638      CALL dom_vvl_interpol( e3v_n(:,:,:), e3vw_n(:,:,:), 'VW' )
639      CALL dom_vvl_interpol( e3t_b(:,:,:),  e3w_b(:,:,:), 'W'  )
640      CALL dom_vvl_interpol( e3u_b(:,:,:), e3uw_b(:,:,:), 'UW' )
641      CALL dom_vvl_interpol( e3v_b(:,:,:), e3vw_b(:,:,:), 'VW' )
[5120]642
[6140]643      ! t- and w- points depth (set the isf depth as it is in the initial step)
[7753]644      gdept_n(:,:,1) = 0.5_wp * e3w_n(:,:,1)
645      gdepw_n(:,:,1) = 0.0_wp
646      gde3w_n(:,:,1) = gdept_n(:,:,1) - sshn(:,:)
[5120]647      DO jk = 2, jpk
648         DO jj = 1,jpj
649            DO ji = 1,jpi
650              !    zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))   ! 0 everywhere tmask = wmask, ie everywhere expect at jk = mikt
651                                                                 ! 1 for jk = mikt
652               zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))
[6140]653               gdepw_n(ji,jj,jk) = gdepw_n(ji,jj,jk-1) + e3t_n(ji,jj,jk-1)
654               gdept_n(ji,jj,jk) =    zcoef  * ( gdepw_n(ji,jj,jk  ) + 0.5 * e3w_n(ji,jj,jk) )  &
655                   &             + (1-zcoef) * ( gdept_n(ji,jj,jk-1) +       e3w_n(ji,jj,jk) ) 
656               gde3w_n(ji,jj,jk) = gdept_n(ji,jj,jk) - sshn(ji,jj)
[4990]657            END DO
658         END DO
[4292]659      END DO
[5120]660
[6140]661      ! Local depth and Inverse of the local depth of the water
662      ! -------------------------------------------------------
[7753]663      hu_n(:,:) = hu_a(:,:)   ;   r1_hu_n(:,:) = r1_hu_a(:,:)
664      hv_n(:,:) = hv_a(:,:)   ;   r1_hv_n(:,:) = r1_hv_a(:,:)
665      !
666      ht_n(:,:) = e3t_n(:,:,1) * tmask(:,:,1)
[6140]667      DO jk = 2, jpkm1
[7753]668         ht_n(:,:) = ht_n(:,:) + e3t_n(:,:,jk) * tmask(:,:,jk)
[4370]669      END DO
[7753]670
[4292]671      ! write restart file
672      ! ==================
[9019]673      IF( lrst_oce  )   CALL dom_vvl_rst( kt, 'WRITE' )
[4292]674      !
[9019]675      IF( ln_timing )   CALL timing_stop('dom_vvl_sf_swp')
[6140]676      !
[4292]677   END SUBROUTINE dom_vvl_sf_swp
678
679
680   SUBROUTINE dom_vvl_interpol( pe3_in, pe3_out, pout )
681      !!---------------------------------------------------------------------
682      !!                  ***  ROUTINE dom_vvl__interpol  ***
683      !!
684      !! ** Purpose :   interpolate scale factors from one grid point to another
685      !!
686      !! ** Method  :   e3_out = e3_0 + interpolation(e3_in - e3_0)
687      !!                - horizontal interpolation: grid cell surface averaging
688      !!                - vertical interpolation: simple averaging
689      !!----------------------------------------------------------------------
[5836]690      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::  pe3_in    ! input e3 to be interpolated
691      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::  pe3_out   ! output interpolated e3
692      CHARACTER(LEN=*)                , INTENT(in   ) ::  pout      ! grid point of out scale factors
693      !                                                             !   =  'U', 'V', 'W, 'F', 'UW' or 'VW'
694      !
[6152]695      INTEGER ::   ji, jj, jk                                       ! dummy loop indices
[9023]696      REAL(wp) ::  zlnwd                                            ! =1./0. when ln_wd_il = T/F
[4292]697      !!----------------------------------------------------------------------
[5836]698      !
[9023]699      IF(ln_wd_il) THEN
[6152]700        zlnwd = 1.0_wp
701      ELSE
702        zlnwd = 0.0_wp
703      END IF
704      !
[5836]705      SELECT CASE ( pout )    !==  type of interpolation  ==!
[4292]706         !
[5836]707      CASE( 'U' )                   !* from T- to U-point : hor. surface weighted mean
[4292]708         DO jk = 1, jpk
709            DO jj = 1, jpjm1
710               DO ji = 1, fs_jpim1   ! vector opt.
[6152]711                  pe3_out(ji,jj,jk) = 0.5_wp * (  umask(ji,jj,jk) * (1.0_wp - zlnwd) + zlnwd ) * r1_e1e2u(ji,jj)   &
[5836]712                     &                       * (   e1e2t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) )     &
713                     &                           + e1e2t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) )
[4292]714               END DO
[2528]715            END DO
716         END DO
[10425]717         CALL lbc_lnk( 'domvvl', pe3_out(:,:,:), 'U', 1._wp )
[7753]718         pe3_out(:,:,:) = pe3_out(:,:,:) + e3u_0(:,:,:)
[5836]719         !
720      CASE( 'V' )                   !* from T- to V-point : hor. surface weighted mean
[4292]721         DO jk = 1, jpk
722            DO jj = 1, jpjm1
723               DO ji = 1, fs_jpim1   ! vector opt.
[6152]724                  pe3_out(ji,jj,jk) = 0.5_wp * ( vmask(ji,jj,jk)  * (1.0_wp - zlnwd) + zlnwd ) * r1_e1e2v(ji,jj)   &
[5836]725                     &                       * (   e1e2t(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3t_0(ji,jj  ,jk) )     &
726                     &                           + e1e2t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) )
[4292]727               END DO
728            END DO
729         END DO
[10425]730         CALL lbc_lnk( 'domvvl', pe3_out(:,:,:), 'V', 1._wp )
[7753]731         pe3_out(:,:,:) = pe3_out(:,:,:) + e3v_0(:,:,:)
[5836]732         !
733      CASE( 'F' )                   !* from U-point to F-point : hor. surface weighted mean
[4292]734         DO jk = 1, jpk
735            DO jj = 1, jpjm1
736               DO ji = 1, fs_jpim1   ! vector opt.
[6152]737                  pe3_out(ji,jj,jk) = 0.5_wp * (  umask(ji,jj,jk) * umask(ji,jj+1,jk) * (1.0_wp - zlnwd) + zlnwd ) &
738                     &                       *    r1_e1e2f(ji,jj)                                                  &
[5836]739                     &                       * (   e1e2u(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3u_0(ji,jj  ,jk) )     &
740                     &                           + e1e2u(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3u_0(ji,jj+1,jk) ) )
[4292]741               END DO
742            END DO
743         END DO
[10425]744         CALL lbc_lnk( 'domvvl', pe3_out(:,:,:), 'F', 1._wp )
[7753]745         pe3_out(:,:,:) = pe3_out(:,:,:) + e3f_0(:,:,:)
[5836]746         !
747      CASE( 'W' )                   !* from T- to W-point : vertical simple mean
748         !
[7753]749         pe3_out(:,:,1) = e3w_0(:,:,1) + pe3_in(:,:,1) - e3t_0(:,:,1)
[5836]750         ! - ML - The use of mask in this formulea enables the special treatment of the last w-point without indirect adressing
751!!gm BUG? use here wmask in case of ISF ?  to be checked
[4292]752         DO jk = 2, jpk
[7753]753            pe3_out(:,:,jk) = e3w_0(:,:,jk) + ( 1.0_wp - 0.5_wp * ( tmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) )   &
754               &                            * ( pe3_in(:,:,jk-1) - e3t_0(:,:,jk-1) )                               &
755               &                            +            0.5_wp * ( tmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd )     &
756               &                            * ( pe3_in(:,:,jk  ) - e3t_0(:,:,jk  ) )
[4292]757         END DO
[5836]758         !
759      CASE( 'UW' )                  !* from U- to UW-point : vertical simple mean
760         !
[7753]761         pe3_out(:,:,1) = e3uw_0(:,:,1) + pe3_in(:,:,1) - e3u_0(:,:,1)
[4292]762         ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without indirect adressing
[5836]763!!gm BUG? use here wumask in case of ISF ?  to be checked
[4292]764         DO jk = 2, jpk
[7753]765            pe3_out(:,:,jk) = e3uw_0(:,:,jk) + ( 1.0_wp - 0.5_wp * ( umask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) )  &
766               &                             * ( pe3_in(:,:,jk-1) - e3u_0(:,:,jk-1) )                              &
767               &                             +            0.5_wp * ( umask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd )    &
768               &                             * ( pe3_in(:,:,jk  ) - e3u_0(:,:,jk  ) )
[4292]769         END DO
[5836]770         !
771      CASE( 'VW' )                  !* from V- to VW-point : vertical simple mean
772         !
[7753]773         pe3_out(:,:,1) = e3vw_0(:,:,1) + pe3_in(:,:,1) - e3v_0(:,:,1)
[4292]774         ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without indirect adressing
[5836]775!!gm BUG? use here wvmask in case of ISF ?  to be checked
[4292]776         DO jk = 2, jpk
[7753]777            pe3_out(:,:,jk) = e3vw_0(:,:,jk) + ( 1.0_wp - 0.5_wp * ( vmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) )  &
778               &                             * ( pe3_in(:,:,jk-1) - e3v_0(:,:,jk-1) )                              &
779               &                             +            0.5_wp * ( vmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd )    &
780               &                             * ( pe3_in(:,:,jk  ) - e3v_0(:,:,jk  ) )
[4292]781         END DO
782      END SELECT
783      !
784   END SUBROUTINE dom_vvl_interpol
785
[5836]786
[4292]787   SUBROUTINE dom_vvl_rst( kt, cdrw )
788      !!---------------------------------------------------------------------
789      !!                   ***  ROUTINE dom_vvl_rst  ***
790      !!                     
791      !! ** Purpose :   Read or write VVL file in restart file
792      !!
793      !! ** Method  :   use of IOM library
794      !!                if the restart does not contain vertical scale factors,
795      !!                they are set to the _0 values
796      !!                if the restart does not contain vertical scale factors increments (z_tilde),
797      !!                they are set to 0.
798      !!----------------------------------------------------------------------
799      INTEGER         , INTENT(in) ::   kt     ! ocean time-step
800      CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
[5836]801      !
[6152]802      INTEGER ::   ji, jj, jk
[4292]803      INTEGER ::   id1, id2, id3, id4, id5     ! local integers
804      !!----------------------------------------------------------------------
805      !
806      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
807         !                                   ! ===============
808         IF( ln_rstart ) THEN                   !* Read the restart file
809            CALL rst_read_open                  !  open the restart file if necessary
[9367]810            CALL iom_get( numror, jpdom_autoglo, 'sshn'   , sshn, ldxios = lrxios    )
[4366]811            !
[6140]812            id1 = iom_varid( numror, 'e3t_b', ldstop = .FALSE. )
813            id2 = iom_varid( numror, 'e3t_n', ldstop = .FALSE. )
[4292]814            id3 = iom_varid( numror, 'tilde_e3t_b', ldstop = .FALSE. )
815            id4 = iom_varid( numror, 'tilde_e3t_n', ldstop = .FALSE. )
[4795]816            id5 = iom_varid( numror, 'hdiv_lf', ldstop = .FALSE. )
[11423]817            !
818            ! coupling with an ice sheet model (grounding/calving line migration)
[11521]819            IF ( ln_isfcpl ) THEN
[11423]820               !
821               ! extrapolate ssh over new water column if needed
822               CALL isfcpl_ssh()
823               !
[11553]824               ! e3t need to be computed where ice shelf draft changed (new ht_0)
[11423]825               id1 = 0 ; id2 = 0 ; id3 = 0 ; id4 = 0 ; id5 = 0
826               !
827            END IF
[4292]828            !                             ! --------- !
829            !                             ! all cases !
830            !                             ! --------- !
831            IF( MIN( id1, id2 ) > 0 ) THEN       ! all required arrays exist
[9367]832               CALL iom_get( numror, jpdom_autoglo, 'e3t_b', e3t_b(:,:,:), ldxios = lrxios )
833               CALL iom_get( numror, jpdom_autoglo, 'e3t_n', e3t_n(:,:,:), ldxios = lrxios )
[4990]834               ! needed to restart if land processor not computed
[6140]835               IF(lwp) write(numout,*) 'dom_vvl_rst : e3t_b and e3t_n found in restart files'
[4990]836               WHERE ( tmask(:,:,:) == 0.0_wp ) 
[6140]837                  e3t_n(:,:,:) = e3t_0(:,:,:)
838                  e3t_b(:,:,:) = e3t_0(:,:,:)
[4990]839               END WHERE
[4292]840               IF( neuler == 0 ) THEN
[6140]841                  e3t_b(:,:,:) = e3t_n(:,:,:)
[4292]842               ENDIF
843            ELSE IF( id1 > 0 ) THEN
[6140]844               IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : e3t_n not found in restart files'
845               IF(lwp) write(numout,*) 'e3t_n set equal to e3t_b.'
[4990]846               IF(lwp) write(numout,*) 'neuler is forced to 0'
[9367]847               CALL iom_get( numror, jpdom_autoglo, 'e3t_b', e3t_b(:,:,:), ldxios = lrxios )
[6140]848               e3t_n(:,:,:) = e3t_b(:,:,:)
[4990]849               neuler = 0
850            ELSE IF( id2 > 0 ) THEN
[6140]851               IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : e3t_b not found in restart files'
852               IF(lwp) write(numout,*) 'e3t_b set equal to e3t_n.'
[4490]853               IF(lwp) write(numout,*) 'neuler is forced to 0'
[9367]854               CALL iom_get( numror, jpdom_autoglo, 'e3t_n', e3t_n(:,:,:), ldxios = lrxios )
[6140]855               e3t_b(:,:,:) = e3t_n(:,:,:)
[4490]856               neuler = 0
857            ELSE
[6140]858               IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : e3t_n not found in restart file'
[4490]859               IF(lwp) write(numout,*) 'Compute scale factor from sshn'
860               IF(lwp) write(numout,*) 'neuler is forced to 0'
[6140]861               DO jk = 1, jpk
862                  e3t_n(:,:,jk) =  e3t_0(:,:,jk) * ( ht_0(:,:) + sshn(:,:) ) &
863                      &                          / ( ht_0(:,:) + 1._wp - ssmask(:,:) ) * tmask(:,:,jk)   &
864                      &          + e3t_0(:,:,jk)                               * (1._wp -tmask(:,:,jk))
[4490]865               END DO
[6140]866               e3t_b(:,:,:) = e3t_n(:,:,:)
[4490]867               neuler = 0
[4292]868            ENDIF
869            !                             ! ----------- !
870            IF( ln_vvl_zstar ) THEN       ! z_star case !
871               !                          ! ----------- !
872               IF( MIN( id3, id4 ) > 0 ) THEN
873                  CALL ctl_stop( 'dom_vvl_rst: z_star cannot restart from a z_tilde or layer run' )
874               ENDIF
875               !                          ! ----------------------- !
876            ELSE                          ! z_tilde and layer cases !
877               !                          ! ----------------------- !
878               IF( MIN( id3, id4 ) > 0 ) THEN  ! all required arrays exist
[9367]879                  CALL iom_get( numror, jpdom_autoglo, 'tilde_e3t_b', tilde_e3t_b(:,:,:), ldxios = lrxios )
880                  CALL iom_get( numror, jpdom_autoglo, 'tilde_e3t_n', tilde_e3t_n(:,:,:), ldxios = lrxios )
[4292]881               ELSE                            ! one at least array is missing
882                  tilde_e3t_b(:,:,:) = 0.0_wp
883                  tilde_e3t_n(:,:,:) = 0.0_wp
884               ENDIF
885               !                          ! ------------ !
886               IF( ln_vvl_ztilde ) THEN   ! z_tilde case !
887                  !                       ! ------------ !
888                  IF( id5 > 0 ) THEN  ! required array exists
[9367]889                     CALL iom_get( numror, jpdom_autoglo, 'hdiv_lf', hdiv_lf(:,:,:), ldxios = lrxios )
[4292]890                  ELSE                ! array is missing
891                     hdiv_lf(:,:,:) = 0.0_wp
892                  ENDIF
893               ENDIF
894            ENDIF
895            !
896         ELSE                                   !* Initialize at "rest"
[7646]897            !
[9023]898
899            IF( ll_wd ) THEN   ! MJB ll_wd edits start here - these are essential
900               !
901               IF( cn_cfg == 'wad' ) THEN
902                  ! Wetting and drying test case
903                  CALL usr_def_istate( gdept_b, tmask, tsb, ub, vb, sshb  )
904                  tsn  (:,:,:,:) = tsb (:,:,:,:)       ! set now values from to before ones
905                  sshn (:,:)     = sshb(:,:)
906                  un   (:,:,:)   = ub  (:,:,:)
907                  vn   (:,:,:)   = vb  (:,:,:)
908               ELSE
909                  ! if not test case
910                  sshn(:,:) = -ssh_ref
911                  sshb(:,:) = -ssh_ref
912
913                  DO jj = 1, jpj
914                     DO ji = 1, jpi
915                        IF( ht_0(ji,jj)-ssh_ref <  rn_wdmin1 ) THEN ! if total depth is less than min depth
916
917                           sshb(ji,jj) = rn_wdmin1 - (ht_0(ji,jj) )
918                           sshn(ji,jj) = rn_wdmin1 - (ht_0(ji,jj) )
919                           ssha(ji,jj) = rn_wdmin1 - (ht_0(ji,jj) )
920                        ENDIF
921                     ENDDO
922                  ENDDO
923               ENDIF !If test case else
924
925               ! Adjust vertical metrics for all wad
[7646]926               DO jk = 1, jpk
[9023]927                  e3t_n(:,:,jk) =  e3t_0(:,:,jk) * ( ht_0(:,:) + sshn(:,:)  ) &
[7646]928                    &                            / ( ht_0(:,:) + 1._wp - ssmask(:,:) ) * tmask(:,:,jk)   &
[9023]929                    &            + e3t_0(:,:,jk) * ( 1._wp - tmask(:,:,jk) )
[7646]930               END DO
931               e3t_b(:,:,:) = e3t_n(:,:,:)
[9023]932
933               DO ji = 1, jpi
934                  DO jj = 1, jpj
935                     IF ( ht_0(ji,jj) .LE. 0.0 .AND. NINT( ssmask(ji,jj) ) .EQ. 1) THEN
936                       CALL ctl_stop( 'dom_vvl_rst: ht_0 must be positive at potentially wet points' )
937                     ENDIF
938                  END DO
939               END DO 
[7646]940               !
941            ELSE
942               !
[9059]943               ! Just to read set ssh in fact, called latter once vertical grid
944               ! is set up:
[9065]945!               CALL usr_def_istate( gdept_0, tmask, tsb, ub, vb, sshb  )
946!               !
947!               DO jk=1,jpk
948!                  e3t_b(:,:,jk) =  e3t_0(:,:,jk) * ( ht_0(:,:) + sshb(:,:) ) &
949!                     &            / ( ht_0(:,:) + 1._wp -ssmask(:,:) ) * tmask(:,:,jk)
950!               END DO
951!               e3t_n(:,:,:) = e3t_b(:,:,:)
952                sshn(:,:)=0._wp
953                e3t_n(:,:,:)=e3t_0(:,:,:)
954                e3t_b(:,:,:)=e3t_0(:,:,:)
[7646]955               !
[9023]956            END IF           ! end of ll_wd edits
[6152]957
[4292]958            IF( ln_vvl_ztilde .OR. ln_vvl_layer) THEN
[7646]959               tilde_e3t_b(:,:,:) = 0._wp
960               tilde_e3t_n(:,:,:) = 0._wp
961               IF( ln_vvl_ztilde ) hdiv_lf(:,:,:) = 0._wp
[4292]962            END IF
963         ENDIF
[5836]964         !
[4292]965      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
966         !                                   ! ===================
967         IF(lwp) WRITE(numout,*) '---- dom_vvl_rst ----'
[9367]968         IF( lwxios ) CALL iom_swap(      cwxios_context          )
[4292]969         !                                           ! --------- !
970         !                                           ! all cases !
971         !                                           ! --------- !
[9367]972         CALL iom_rstput( kt, nitrst, numrow, 'e3t_b', e3t_b(:,:,:), ldxios = lwxios )
973         CALL iom_rstput( kt, nitrst, numrow, 'e3t_n', e3t_n(:,:,:), ldxios = lwxios )
[4292]974         !                                           ! ----------------------- !
975         IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN  ! z_tilde and layer cases !
976            !                                        ! ----------------------- !
[9367]977            CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_b', tilde_e3t_b(:,:,:), ldxios = lwxios)
978            CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_n', tilde_e3t_n(:,:,:), ldxios = lwxios)
[4292]979         END IF
980         !                                           ! -------------!   
981         IF( ln_vvl_ztilde ) THEN                    ! z_tilde case !
982            !                                        ! ------------ !
[9367]983            CALL iom_rstput( kt, nitrst, numrow, 'hdiv_lf', hdiv_lf(:,:,:), ldxios = lwxios)
[4292]984         ENDIF
[5836]985         !
[9367]986         IF( lwxios ) CALL iom_swap(      cxios_context          )
[4292]987      ENDIF
[5836]988      !
[4292]989   END SUBROUTINE dom_vvl_rst
990
991
992   SUBROUTINE dom_vvl_ctl
993      !!---------------------------------------------------------------------
994      !!                  ***  ROUTINE dom_vvl_ctl  ***
995      !!               
996      !! ** Purpose :   Control the consistency between namelist options
997      !!                for vertical coordinate
998      !!----------------------------------------------------------------------
[5836]999      INTEGER ::   ioptio, ios
1000      !!
[4292]1001      NAMELIST/nam_vvl/ ln_vvl_zstar, ln_vvl_ztilde, ln_vvl_layer, ln_vvl_ztilde_as_zstar, &
[5836]1002         &              ln_vvl_zstar_at_eqtor      , rn_ahe3     , rn_rst_e3t            , &
1003         &              rn_lf_cutoff               , rn_zdef_max , ln_vvl_dbg                ! not yet implemented: ln_vvl_kepe
[4292]1004      !!----------------------------------------------------------------------
[5836]1005      !
[4294]1006      REWIND( numnam_ref )              ! Namelist nam_vvl in reference namelist :
1007      READ  ( numnam_ref, nam_vvl, IOSTAT = ios, ERR = 901)
[9168]1008901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nam_vvl in reference namelist', lwp )
[4294]1009      REWIND( numnam_cfg )              ! Namelist nam_vvl in configuration namelist : Parameters of the run
1010      READ  ( numnam_cfg, nam_vvl, IOSTAT = ios, ERR = 902 )
[9168]1011902   IF( ios >  0 ) CALL ctl_nam ( ios , 'nam_vvl in configuration namelist', lwp )
[4624]1012      IF(lwm) WRITE ( numond, nam_vvl )
[5836]1013      !
[4292]1014      IF(lwp) THEN                    ! Namelist print
1015         WRITE(numout,*)
1016         WRITE(numout,*) 'dom_vvl_ctl : choice/control of the variable vertical coordinate'
1017         WRITE(numout,*) '~~~~~~~~~~~'
[9168]1018         WRITE(numout,*) '   Namelist nam_vvl : chose a vertical coordinate'
1019         WRITE(numout,*) '      zstar                      ln_vvl_zstar   = ', ln_vvl_zstar
1020         WRITE(numout,*) '      ztilde                     ln_vvl_ztilde  = ', ln_vvl_ztilde
1021         WRITE(numout,*) '      layer                      ln_vvl_layer   = ', ln_vvl_layer
1022         WRITE(numout,*) '      ztilde as zstar   ln_vvl_ztilde_as_zstar  = ', ln_vvl_ztilde_as_zstar
[4292]1023         WRITE(numout,*) '      ztilde near the equator    ln_vvl_zstar_at_eqtor  = ', ln_vvl_zstar_at_eqtor
[9168]1024         WRITE(numout,*) '      !'
1025         WRITE(numout,*) '      thickness diffusion coefficient                      rn_ahe3      = ', rn_ahe3
1026         WRITE(numout,*) '      maximum e3t deformation fractional change            rn_zdef_max  = ', rn_zdef_max
[4292]1027         IF( ln_vvl_ztilde_as_zstar ) THEN
[9168]1028            WRITE(numout,*) '      ztilde running in zstar emulation mode (ln_vvl_ztilde_as_zstar=T) '
1029            WRITE(numout,*) '         ignoring namelist timescale parameters and using:'
1030            WRITE(numout,*) '            hard-wired : z-tilde to zstar restoration timescale (days)'
1031            WRITE(numout,*) '                         rn_rst_e3t     = 0.e0'
1032            WRITE(numout,*) '            hard-wired : z-tilde cutoff frequency of low-pass filter (days)'
1033            WRITE(numout,*) '                         rn_lf_cutoff   = 1.0/rdt'
[4292]1034         ELSE
[9168]1035            WRITE(numout,*) '      z-tilde to zstar restoration timescale (days)        rn_rst_e3t   = ', rn_rst_e3t
1036            WRITE(numout,*) '      z-tilde cutoff frequency of low-pass filter (days)   rn_lf_cutoff = ', rn_lf_cutoff
[4292]1037         ENDIF
[9168]1038         WRITE(numout,*) '         debug prints flag                                 ln_vvl_dbg   = ', ln_vvl_dbg
[4292]1039      ENDIF
[5836]1040      !
[4292]1041      ioptio = 0                      ! Parameter control
[5836]1042      IF( ln_vvl_ztilde_as_zstar )   ln_vvl_ztilde = .true.
1043      IF( ln_vvl_zstar           )   ioptio = ioptio + 1
1044      IF( ln_vvl_ztilde          )   ioptio = ioptio + 1
1045      IF( ln_vvl_layer           )   ioptio = ioptio + 1
1046      !
[4292]1047      IF( ioptio /= 1 )   CALL ctl_stop( 'Choose ONE vertical coordinate in namelist nam_vvl' )
[6140]1048      IF( .NOT. ln_vvl_zstar .AND. ln_isf ) CALL ctl_stop( 'Only vvl_zstar has been tested with ice shelf cavity' )
[5836]1049      !
[4292]1050      IF(lwp) THEN                   ! Print the choice
1051         WRITE(numout,*)
[9168]1052         IF( ln_vvl_zstar           ) WRITE(numout,*) '      ==>>>   zstar vertical coordinate is used'
1053         IF( ln_vvl_ztilde          ) WRITE(numout,*) '      ==>>>   ztilde vertical coordinate is used'
1054         IF( ln_vvl_layer           ) WRITE(numout,*) '      ==>>>   layer vertical coordinate is used'
1055         IF( ln_vvl_ztilde_as_zstar ) WRITE(numout,*) '      ==>>>   to emulate a zstar coordinate'
[4292]1056      ENDIF
[5836]1057      !
[4486]1058#if defined key_agrif
[9190]1059      IF( (.NOT.Agrif_Root()).AND.(.NOT.ln_vvl_zstar) )   CALL ctl_stop( 'AGRIF is implemented with zstar coordinate only' )
[4486]1060#endif
[5836]1061      !
[4292]1062   END SUBROUTINE dom_vvl_ctl
1063
[592]1064   !!======================================================================
1065END MODULE domvvl
Note: See TracBrowser for help on using the repository browser.