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/dev_r11943_MERGE_2019/src/OCE/DOM – NEMO

source: NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/DOM/domvvl.F90 @ 11949

Last change on this file since 11949 was 11949, checked in by acc, 5 years ago

Merge in changes from 2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps. This just creates a fresh copy of this branch to use as the merge base. See ticket #2341

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