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

source: NEMO/trunk/src/OCE/DOM/domvvl.F90

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

#1791: correct thicknesses at last level with vvl. This ensures total depth is correct (eg matches depth at rest + sea level anomaly) which may otherwise induce unstabilities, and unconsistencies with external mode. The correction is active over steps (no change with purely s-coordinates or in between two adjacent z-levels with partial cells). Contrary to qco formulation, the vvl scheme can lead to negative thicknesses at velocity points, even if total depth at T-points is positive. This is a rare case which is not taken into account in this correction.

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