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_r13723_KERNEL-01_Amy_Mike_newHPGschemes/src/OCE/DOM – NEMO

source: NEMO/branches/2020/dev_r13723_KERNEL-01_Amy_Mike_newHPGschemes/src/OCE/DOM/domvvl.F90 @ 15475

Last change on this file since 15475 was 14058, checked in by ayoung, 4 years ago

Updating to trunk revision 14057. No conflicts since last sette test. Ticket #2480.

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