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/2021/dev_r14318_RK3_stage1_tsplit/tests/SEICHE/MY_SRC – NEMO

source: NEMO/branches/2021/dev_r14318_RK3_stage1_tsplit/tests/SEICHE/MY_SRC/domvvl.F90 @ 14597

Last change on this file since 14597 was 14597, checked in by jchanut, 3 years ago

#2634: Implement Marsaleix (2008) and Demange (2019) test case (an external SEICHE over a 2DV stratified ocean and a seamount). One can optionnaly retrieve a standard external or 2 layer internal Seiche with a flat bottom.

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