New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
domvvl.F90 in NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/DOM – NEMO

source: NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/src/OCE/DOM/domvvl.F90 @ 10743

Last change on this file since 10743 was 10743, checked in by davestorkey, 5 years ago

First block of changes:

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