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/2020/dev_r13312_AGRIF-03-04_jchanut_vinterp_tstep/src/OCE/DOM – NEMO

source: NEMO/branches/2020/dev_r13312_AGRIF-03-04_jchanut_vinterp_tstep/src/OCE/DOM/domvvl.F90 @ 13334

Last change on this file since 13334 was 13334, checked in by jchanut, 4 years ago

finish bypassing ocean/ice initialization with AGRIF, #2222, #2129

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