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 branches/UKMO/dev_r5518_GO6_package_XIOS_read/NEMOGCM/NEMO/OPA_SRC/DOM – NEMO

source: branches/UKMO/dev_r5518_GO6_package_XIOS_read/NEMOGCM/NEMO/OPA_SRC/DOM/domvvl.F90 @ 7924

Last change on this file since 7924 was 7924, checked in by andmirek, 7 years ago

first commit with XIOS restart read functionality

File size: 75.3 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:
9   !!                                          vvl option includes z_star and z_tilde coordinates
10   !!            3.6  !  2014-11  (P. Mathiot) add ice shelf capability
11   !!----------------------------------------------------------------------
12   !!   'key_vvl'                              variable volume
13   !!----------------------------------------------------------------------
14   !!----------------------------------------------------------------------
15   !!   dom_vvl_init     : define initial vertical scale factors, depths and column thickness
16   !!   dom_vvl_sf_nxt   : Compute next vertical scale factors
17   !!   dom_vvl_sf_swp   : Swap vertical scale factors and update the vertical grid
18   !!   dom_vvl_interpol : Interpolate vertical scale factors from one grid point to another
19   !!   dom_vvl_rst      : read/write restart file
20   !!   dom_vvl_ctl      : Check the vvl options
21   !!   dom_vvl_orca_fix : Recompute some area-weighted interpolations of vertical scale factors
22   !!                    : to account for manual changes to e[1,2][u,v] in some Straits
23   !!----------------------------------------------------------------------
24   !! * Modules used
25   USE oce             ! ocean dynamics and tracers
26   USE dom_oce         ! ocean space and time domain
27   USE sbc_oce         ! ocean surface boundary condition
28   USE in_out_manager  ! I/O manager
29   USE iom             ! I/O manager library
30   USE restart         ! ocean restart
31   USE lib_mpp         ! distributed memory computing library
32   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
33   USE wrk_nemo        ! Memory allocation
34   USE timing          ! Timing
35   USE iom_def, ONLY : lxios_read
36
37   IMPLICIT NONE
38   PRIVATE
39
40   !! * Routine accessibility
41   PUBLIC  dom_vvl_init       ! called by domain.F90
42   PUBLIC  dom_vvl_sf_nxt     ! called by step.F90
43   PUBLIC  dom_vvl_sf_swp     ! called by step.F90
44   PUBLIC  dom_vvl_interpol   ! called by dynnxt.F90
45   PRIVATE dom_vvl_orca_fix   ! called by dom_vvl_interpol
46
47   !!* Namelist nam_vvl
48   LOGICAL , PUBLIC                                      :: ln_vvl_zstar = .FALSE.              ! zstar  vertical coordinate
49   LOGICAL , PUBLIC                                      :: ln_vvl_ztilde = .FALSE.             ! ztilde vertical coordinate
50   LOGICAL , PUBLIC                                      :: ln_vvl_layer = .FALSE.              ! level  vertical coordinate
51   LOGICAL , PUBLIC                                      :: ln_vvl_ztilde_as_zstar = .FALSE.    ! ztilde vertical coordinate
52   LOGICAL , PUBLIC                                      :: ln_vvl_zstar_at_eqtor = .FALSE.     ! ztilde vertical coordinate
53   LOGICAL , PUBLIC                                      :: ln_vvl_kepe = .FALSE.               ! kinetic/potential energy transfer
54   !                                                                                            ! conservation: not used yet
55   REAL(wp)                                              :: rn_ahe3                   ! thickness diffusion coefficient
56   REAL(wp)                                              :: rn_rst_e3t                ! ztilde to zstar restoration timescale [days]
57   REAL(wp)                                              :: rn_lf_cutoff              ! cutoff frequency for low-pass filter  [days]
58   REAL(wp)                                              :: rn_zdef_max               ! maximum fractional e3t deformation
59   LOGICAL , PUBLIC                                      :: ln_vvl_dbg = .FALSE.      ! debug control prints
60
61   !! * Module variables
62   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: un_td, vn_td                       ! thickness diffusion transport
63   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: hdiv_lf                            ! low frequency part of hz divergence
64   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: tilde_e3t_b, tilde_e3t_n           ! baroclinic scale factors
65   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: tilde_e3t_a, dtilde_e3t_a          ! baroclinic scale factors
66   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   :: frq_rst_e3t                        ! retoring period for scale factors
67   REAL(wp)        , ALLOCATABLE, SAVE, DIMENSION(:,:)   :: frq_rst_hdv                        ! retoring period for low freq. divergence
68
69   !! * Substitutions
70#  include "domzgr_substitute.h90"
71#  include "vectopt_loop_substitute.h90"
72   !!----------------------------------------------------------------------
73   !! NEMO/OPA 3.3 , NEMO-Consortium (2010)
74   !! $Id$
75   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
76   !!----------------------------------------------------------------------
77
78CONTAINS
79
80   INTEGER FUNCTION dom_vvl_alloc()
81      !!----------------------------------------------------------------------
82      !!                ***  FUNCTION dom_vvl_alloc  ***
83      !!----------------------------------------------------------------------
84      IF( ln_vvl_zstar ) dom_vvl_alloc = 0
85      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
86         ALLOCATE( tilde_e3t_b(jpi,jpj,jpk)  , tilde_e3t_n(jpi,jpj,jpk) , tilde_e3t_a(jpi,jpj,jpk) ,   &
87            &      dtilde_e3t_a(jpi,jpj,jpk) , un_td  (jpi,jpj,jpk)     , vn_td  (jpi,jpj,jpk)     ,   &
88            &      STAT = dom_vvl_alloc        )
89         IF( lk_mpp             )   CALL mpp_sum ( dom_vvl_alloc )
90         IF( dom_vvl_alloc /= 0 )   CALL ctl_warn('dom_vvl_alloc: failed to allocate arrays')
91         un_td = 0.0_wp
92         vn_td = 0.0_wp
93      ENDIF
94      IF( ln_vvl_ztilde ) THEN
95         ALLOCATE( frq_rst_e3t(jpi,jpj) , frq_rst_hdv(jpi,jpj) , hdiv_lf(jpi,jpj,jpk) , STAT= dom_vvl_alloc )
96         IF( lk_mpp             )   CALL mpp_sum ( dom_vvl_alloc )
97         IF( dom_vvl_alloc /= 0 )   CALL ctl_warn('dom_vvl_alloc: failed to allocate arrays')
98      ENDIF
99
100   END FUNCTION dom_vvl_alloc
101
102
103   SUBROUTINE dom_vvl_init
104      !!----------------------------------------------------------------------
105      !!                ***  ROUTINE dom_vvl_init  ***
106      !!                   
107      !! ** Purpose :  Initialization of all scale factors, depths
108      !!               and water column heights
109      !!
110      !! ** Method  :  - use restart file and/or initialize
111      !!               - interpolate scale factors
112      !!
113      !! ** Action  : - fse3t_(n/b) and tilde_e3t_(n/b)
114      !!              - Regrid: fse3(u/v)_n
115      !!                        fse3(u/v)_b       
116      !!                        fse3w_n           
117      !!                        fse3(u/v)w_b     
118      !!                        fse3(u/v)w_n     
119      !!                        fsdept_n, fsdepw_n and fsde3w_n
120      !!              - h(t/u/v)_0
121      !!              - frq_rst_e3t and frq_rst_hdv
122      !!
123      !! Reference  : Leclair, M., and G. Madec, 2011, Ocean Modelling.
124      !!----------------------------------------------------------------------
125      USE phycst,  ONLY : rpi, rsmall, rad
126      !! * Local declarations
127      INTEGER ::   ji,jj,jk
128      INTEGER ::   ii0, ii1, ij0, ij1
129      REAL(wp)::   zcoef
130      !!----------------------------------------------------------------------
131      IF( nn_timing == 1 )  CALL timing_start('dom_vvl_init')
132
133      IF(lwp) WRITE(numout,*)
134      IF(lwp) WRITE(numout,*) 'dom_vvl_init : Variable volume activated'
135      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~'
136
137      ! choose vertical coordinate (z_star, z_tilde or layer)
138      ! ==========================
139      CALL dom_vvl_ctl
140
141      ! Allocate module arrays
142      ! ======================
143      IF( dom_vvl_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'dom_vvl_init : unable to allocate arrays' )
144
145      ! Read or initialize fse3t_(b/n), tilde_e3t_(b/n) and hdiv_lf (and e3t_a(jpk))
146      ! ============================================================================
147      CALL dom_vvl_rst( nit000, 'READ' )
148      fse3t_a(:,:,jpk) = e3t_0(:,:,jpk)
149
150      ! Reconstruction of all vertical scale factors at now and before time steps
151      ! =============================================================================
152      ! Horizontal scale factor interpolations
153      ! --------------------------------------
154      CALL dom_vvl_interpol( fse3t_b(:,:,:), fse3u_b(:,:,:), 'U' )
155      CALL dom_vvl_interpol( fse3t_b(:,:,:), fse3v_b(:,:,:), 'V' )
156      CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3u_n(:,:,:), 'U' )
157      CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3v_n(:,:,:), 'V' )
158      CALL dom_vvl_interpol( fse3u_n(:,:,:), fse3f_n(:,:,:), 'F' )
159      ! Vertical scale factor interpolations
160      ! ------------------------------------
161      CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3w_n (:,:,:), 'W'  )
162      CALL dom_vvl_interpol( fse3u_n(:,:,:), fse3uw_n(:,:,:), 'UW' )
163      CALL dom_vvl_interpol( fse3v_n(:,:,:), fse3vw_n(:,:,:), 'VW' )
164      CALL dom_vvl_interpol( fse3t_b(:,:,:), fse3w_b (:,:,:), 'W'  )
165      CALL dom_vvl_interpol( fse3u_b(:,:,:), fse3uw_b(:,:,:), 'UW' )
166      CALL dom_vvl_interpol( fse3v_b(:,:,:), fse3vw_b(:,:,:), 'VW' )
167      ! t- and w- points depth
168      ! ----------------------
169      ! set the isf depth as it is in the initial step
170      fsdept_n(:,:,1) = 0.5_wp * fse3w_n(:,:,1)
171      fsdepw_n(:,:,1) = 0.0_wp
172      fsde3w_n(:,:,1) = fsdept_n(:,:,1) - sshn(:,:)
173      fsdept_b(:,:,1) = 0.5_wp * fse3w_b(:,:,1)
174      fsdepw_b(:,:,1) = 0.0_wp
175
176      DO jk = 2, jpk
177         DO jj = 1,jpj
178            DO ji = 1,jpi
179              !    zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))   ! 0 everywhere tmask = wmask, ie everywhere expect at jk = mikt
180                                                     ! 1 everywhere from mbkt to mikt + 1 or 1 (if no isf)
181                                                     ! 0.5 where jk = mikt 
182               zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))
183               fsdepw_n(ji,jj,jk) = fsdepw_n(ji,jj,jk-1) + fse3t_n(ji,jj,jk-1)
184               fsdept_n(ji,jj,jk) =      zcoef  * ( fsdepw_n(ji,jj,jk  ) + 0.5 * fse3w_n(ji,jj,jk))  &
185                   &                + (1-zcoef) * ( fsdept_n(ji,jj,jk-1) +       fse3w_n(ji,jj,jk)) 
186               fsde3w_n(ji,jj,jk) = fsdept_n(ji,jj,jk) - sshn(ji,jj)
187               fsdepw_b(ji,jj,jk) = fsdepw_b(ji,jj,jk-1) + fse3t_b(ji,jj,jk-1)
188               fsdept_b(ji,jj,jk) =      zcoef  * ( fsdepw_b(ji,jj,jk  ) + 0.5 * fse3w_b(ji,jj,jk))  &
189                   &                + (1-zcoef) * ( fsdept_b(ji,jj,jk-1) +       fse3w_b(ji,jj,jk)) 
190            END DO
191         END DO
192      END DO
193
194      ! Before depth and Inverse of the local depth of the water column at u- and v- points
195      ! -----------------------------------------------------------------------------------
196      hu_b(:,:) = 0.
197      hv_b(:,:) = 0.
198      DO jk = 1, jpkm1
199         hu_b(:,:) = hu_b(:,:) + fse3u_b(:,:,jk) * umask(:,:,jk)
200         hv_b(:,:) = hv_b(:,:) + fse3v_b(:,:,jk) * vmask(:,:,jk)
201      END DO
202      hur_b(:,:) = umask_i(:,:) / ( hu_b(:,:) + 1. - umask_i(:,:) )
203      hvr_b(:,:) = vmask_i(:,:) / ( hv_b(:,:) + 1. - vmask_i(:,:) )
204
205      ! Restoring frequencies for z_tilde coordinate
206      ! ============================================
207      IF( ln_vvl_ztilde ) THEN
208         ! Values in days provided via the namelist; use rsmall to avoid possible division by zero errors with faulty settings
209         frq_rst_e3t(:,:) = 2.0_wp * rpi / ( MAX( rn_rst_e3t  , rsmall ) * 86400.0_wp )
210         frq_rst_hdv(:,:) = 2.0_wp * rpi / ( MAX( rn_lf_cutoff, rsmall ) * 86400.0_wp )
211         IF( ln_vvl_ztilde_as_zstar ) THEN
212            ! Ignore namelist settings and use these next two to emulate z-star using z-tilde
213            frq_rst_e3t(:,:) = 0.0_wp 
214            frq_rst_hdv(:,:) = 1.0_wp / rdt
215         ENDIF
216         IF ( ln_vvl_zstar_at_eqtor ) THEN
217            DO jj = 1, jpj
218               DO ji = 1, jpi
219                  IF( ABS(gphit(ji,jj)) >= 6.) THEN
220                     ! values outside the equatorial band and transition zone (ztilde)
221                     frq_rst_e3t(ji,jj) =  2.0_wp * rpi / ( MAX( rn_rst_e3t  , rsmall ) * 86400.e0_wp )
222                     frq_rst_hdv(ji,jj) =  2.0_wp * rpi / ( MAX( rn_lf_cutoff, rsmall ) * 86400.e0_wp )
223                  ELSEIF( ABS(gphit(ji,jj)) <= 2.5) THEN
224                     ! values inside the equatorial band (ztilde as zstar)
225                     frq_rst_e3t(ji,jj) =  0.0_wp
226                     frq_rst_hdv(ji,jj) =  1.0_wp / rdt
227                  ELSE
228                     ! values in the transition band (linearly vary from ztilde to ztilde as zstar values)
229                     frq_rst_e3t(ji,jj) = 0.0_wp + (frq_rst_e3t(ji,jj)-0.0_wp)*0.5_wp   &
230                        &            * (  1.0_wp - COS( rad*(ABS(gphit(ji,jj))-2.5_wp)  &
231                        &                                          * 180._wp / 3.5_wp ) )
232                     frq_rst_hdv(ji,jj) = (1.0_wp / rdt)                                &
233                        &            + (  frq_rst_hdv(ji,jj)-(1.e0_wp / rdt) )*0.5_wp   &
234                        &            * (  1._wp  - COS( rad*(ABS(gphit(ji,jj))-2.5_wp)  &
235                        &                                          * 180._wp / 3.5_wp ) )
236                  ENDIF
237               END DO
238            END DO
239            IF( cp_cfg == "orca" .AND. jp_cfg == 3 ) THEN
240               ii0 = 103   ;   ii1 = 111        ! Suppress ztilde in the Foxe Basin for ORCA2
241               ij0 = 128   ;   ij1 = 135   ;   
242               frq_rst_e3t( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) =  0.0_wp
243               frq_rst_hdv( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) =  1.e0_wp / rdt
244            ENDIF
245         ENDIF
246      ENDIF
247
248      IF( nn_timing == 1 )  CALL timing_stop('dom_vvl_init')
249
250   END SUBROUTINE dom_vvl_init
251
252
253   SUBROUTINE dom_vvl_sf_nxt( kt, kcall ) 
254      !!----------------------------------------------------------------------
255      !!                ***  ROUTINE dom_vvl_sf_nxt  ***
256      !!                   
257      !! ** Purpose :  - compute the after scale factors used in tra_zdf, dynnxt,
258      !!                 tranxt and dynspg routines
259      !!
260      !! ** Method  :  - z_star case:  Repartition of ssh INCREMENT proportionnaly to the level thickness.
261      !!               - z_tilde_case: after scale factor increment =
262      !!                                    high frequency part of horizontal divergence
263      !!                                  + retsoring towards the background grid
264      !!                                  + thickness difusion
265      !!                               Then repartition of ssh INCREMENT proportionnaly
266      !!                               to the "baroclinic" level thickness.
267      !!
268      !! ** Action  :  - hdiv_lf    : restoring towards full baroclinic divergence in z_tilde case
269      !!               - tilde_e3t_a: after increment of vertical scale factor
270      !!                              in z_tilde case
271      !!               - fse3(t/u/v)_a
272      !!
273      !! Reference  : Leclair, M., and Madec, G. 2011, Ocean Modelling.
274      !!----------------------------------------------------------------------
275      REAL(wp), POINTER, DIMENSION(:,:,:) :: ze3t
276      REAL(wp), POINTER, DIMENSION(:,:  ) :: zht, z_scale, zwu, zwv, zhdiv
277      !! * Arguments
278      INTEGER, INTENT( in )                  :: kt                    ! time step
279      INTEGER, INTENT( in ), OPTIONAL        :: kcall                 ! optional argument indicating call sequence
280      !! * Local declarations
281      INTEGER                                :: ji, jj, jk            ! dummy loop indices
282      INTEGER , DIMENSION(3)                 :: ijk_max, ijk_min      ! temporary integers
283      REAL(wp)                               :: z2dt                  ! temporary scalars
284      REAL(wp)                               :: z_tmin, z_tmax        ! temporary scalars
285      LOGICAL                                :: ll_do_bclinic         ! temporary logical
286      !!----------------------------------------------------------------------
287      IF( nn_timing == 1 )  CALL timing_start('dom_vvl_sf_nxt')
288      CALL wrk_alloc( jpi, jpj, zht, z_scale, zwu, zwv, zhdiv )
289      CALL wrk_alloc( jpi, jpj, jpk, ze3t                     )
290
291      IF(kt == nit000)   THEN
292         IF(lwp) WRITE(numout,*)
293         IF(lwp) WRITE(numout,*) 'dom_vvl_sf_nxt : compute after scale factors'
294         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~'
295      ENDIF
296
297      ll_do_bclinic = .TRUE.
298      IF( PRESENT(kcall) ) THEN
299         IF ( kcall == 2 .AND. ln_vvl_ztilde ) ll_do_bclinic = .FALSE.
300      ENDIF
301
302      ! ******************************* !
303      ! After acale factors at t-points !
304      ! ******************************* !
305
306      !                                                ! --------------------------------------------- !
307                                                       ! z_star coordinate and barotropic z-tilde part !
308      !                                                ! --------------------------------------------- !
309
310      z_scale(:,:) = ( ssha(:,:) - sshb(:,:) ) * ssmask(:,:) / ( ht_0(:,:) + sshn(:,:) + 1. - ssmask(:,:) )
311      DO jk = 1, jpkm1
312         ! formally this is the same as fse3t_a = e3t_0*(1+ssha/ht_0)
313         fse3t_a(:,:,jk) = fse3t_b(:,:,jk) + fse3t_n(:,:,jk) * z_scale(:,:) * tmask(:,:,jk)
314      END DO
315
316      IF( ln_vvl_ztilde .OR. ln_vvl_layer .AND. ll_do_bclinic ) THEN   ! z_tilde or layer coordinate !
317         !                                                            ! ------baroclinic part------ !
318
319         ! I - initialization
320         ! ==================
321
322         ! 1 - barotropic divergence
323         ! -------------------------
324         zhdiv(:,:) = 0.
325         zht(:,:)   = 0.
326         DO jk = 1, jpkm1
327            zhdiv(:,:) = zhdiv(:,:) + fse3t_n(:,:,jk) * hdivn(:,:,jk)
328            zht  (:,:) = zht  (:,:) + fse3t_n(:,:,jk) * tmask(:,:,jk)
329         END DO
330         zhdiv(:,:) = zhdiv(:,:) / ( zht(:,:) + 1. - tmask_i(:,:) )
331
332         ! 2 - Low frequency baroclinic horizontal divergence  (z-tilde case only)
333         ! --------------------------------------------------
334         IF( ln_vvl_ztilde ) THEN
335            IF( kt .GT. nit000 ) THEN
336               DO jk = 1, jpkm1
337                  hdiv_lf(:,:,jk) = hdiv_lf(:,:,jk) - rdt * frq_rst_hdv(:,:)   &
338                     &          * ( hdiv_lf(:,:,jk) - fse3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) )
339               END DO
340            ENDIF
341         END IF
342
343         ! II - after z_tilde increments of vertical scale factors
344         ! =======================================================
345         tilde_e3t_a(:,:,:) = 0.0_wp  ! tilde_e3t_a used to store tendency terms
346
347         ! 1 - High frequency divergence term
348         ! ----------------------------------
349         IF( ln_vvl_ztilde ) THEN     ! z_tilde case
350            DO jk = 1, jpkm1
351               tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - ( fse3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) - hdiv_lf(:,:,jk) )
352            END DO
353         ELSE                         ! layer case
354            DO jk = 1, jpkm1
355               tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) -   fse3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) * tmask(:,:,jk)
356            END DO
357         END IF
358
359         ! 2 - Restoring term (z-tilde case only)
360         ! ------------------
361         IF( ln_vvl_ztilde ) THEN
362            DO jk = 1, jpk
363               tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - frq_rst_e3t(:,:) * tilde_e3t_b(:,:,jk)
364            END DO
365         END IF
366
367         ! 3 - Thickness diffusion term
368         ! ----------------------------
369         zwu(:,:) = 0.0_wp
370         zwv(:,:) = 0.0_wp
371         ! a - first derivative: diffusive fluxes
372         DO jk = 1, jpkm1
373            DO jj = 1, jpjm1
374               DO ji = 1, fs_jpim1   ! vector opt.
375                  un_td(ji,jj,jk) = rn_ahe3 * umask(ji,jj,jk) * re2u_e1u(ji,jj) &
376                                  & * ( tilde_e3t_b(ji,jj,jk) - tilde_e3t_b(ji+1,jj  ,jk) )
377                  vn_td(ji,jj,jk) = rn_ahe3 * vmask(ji,jj,jk) * re1v_e2v(ji,jj) & 
378                                  & * ( tilde_e3t_b(ji,jj,jk) - tilde_e3t_b(ji  ,jj+1,jk) )
379                  zwu(ji,jj) = zwu(ji,jj) + un_td(ji,jj,jk)
380                  zwv(ji,jj) = zwv(ji,jj) + vn_td(ji,jj,jk)
381               END DO
382            END DO
383         END DO
384         ! b - correction for last oceanic u-v points
385         DO jj = 1, jpj
386            DO ji = 1, jpi
387               un_td(ji,jj,mbku(ji,jj)) = un_td(ji,jj,mbku(ji,jj)) - zwu(ji,jj)
388               vn_td(ji,jj,mbkv(ji,jj)) = vn_td(ji,jj,mbkv(ji,jj)) - zwv(ji,jj)
389            END DO
390         END DO
391         ! c - second derivative: divergence of diffusive fluxes
392         DO jk = 1, jpkm1
393            DO jj = 2, jpjm1
394               DO ji = fs_2, fs_jpim1   ! vector opt.
395                  tilde_e3t_a(ji,jj,jk) = tilde_e3t_a(ji,jj,jk) + (   un_td(ji-1,jj  ,jk) - un_td(ji,jj,jk)    &
396                     &                                          +     vn_td(ji  ,jj-1,jk) - vn_td(ji,jj,jk)    &
397                     &                                            ) * r1_e12t(ji,jj)
398               END DO
399            END DO
400         END DO
401         ! d - thickness diffusion transport: boundary conditions
402         !     (stored for tracer advction and continuity equation)
403         CALL lbc_lnk( un_td , 'U' , -1._wp)
404         CALL lbc_lnk( vn_td , 'V' , -1._wp)
405
406         ! 4 - Time stepping of baroclinic scale factors
407         ! ---------------------------------------------
408         ! Leapfrog time stepping
409         ! ~~~~~~~~~~~~~~~~~~~~~~
410         IF( neuler == 0 .AND. kt == nit000 ) THEN
411            z2dt =  rdt
412         ELSE
413            z2dt = 2.0_wp * rdt
414         ENDIF
415         CALL lbc_lnk( tilde_e3t_a(:,:,:), 'T', 1._wp )
416         tilde_e3t_a(:,:,:) = tilde_e3t_b(:,:,:) + z2dt * tmask(:,:,:) * tilde_e3t_a(:,:,:)
417
418         ! Maximum deformation control
419         ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~
420         ze3t(:,:,jpk) = 0.0_wp
421         DO jk = 1, jpkm1
422            ze3t(:,:,jk) = tilde_e3t_a(:,:,jk) / e3t_0(:,:,jk) * tmask(:,:,jk) * tmask_i(:,:)
423         END DO
424         z_tmax = MAXVAL( ze3t(:,:,:) )
425         IF( lk_mpp )   CALL mpp_max( z_tmax )                 ! max over the global domain
426         z_tmin = MINVAL( ze3t(:,:,:) )
427         IF( lk_mpp )   CALL mpp_min( z_tmin )                 ! min over the global domain
428         ! - ML - test: for the moment, stop simulation for too large e3_t variations
429         IF( ( z_tmax .GT. rn_zdef_max ) .OR. ( z_tmin .LT. - rn_zdef_max ) ) THEN
430            IF( lk_mpp ) THEN
431               CALL mpp_maxloc( ze3t, tmask, z_tmax, ijk_max(1), ijk_max(2), ijk_max(3) )
432               CALL mpp_minloc( ze3t, tmask, z_tmin, ijk_min(1), ijk_min(2), ijk_min(3) )
433            ELSE
434               ijk_max = MAXLOC( ze3t(:,:,:) )
435               ijk_max(1) = ijk_max(1) + nimpp - 1
436               ijk_max(2) = ijk_max(2) + njmpp - 1
437               ijk_min = MINLOC( ze3t(:,:,:) )
438               ijk_min(1) = ijk_min(1) + nimpp - 1
439               ijk_min(2) = ijk_min(2) + njmpp - 1
440            ENDIF
441            IF (lwp) THEN
442               WRITE(numout, *) 'MAX( tilde_e3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmax
443               WRITE(numout, *) 'at i, j, k=', ijk_max
444               WRITE(numout, *) 'MIN( tilde_e3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmin
445               WRITE(numout, *) 'at i, j, k=', ijk_min           
446               CALL ctl_warn('MAX( ABS( tilde_e3t_a(:,:,:) ) / e3t_0(:,:,:) ) too high')
447            ENDIF
448         ENDIF
449         ! - ML - end test
450         ! - ML - Imposing these limits will cause a baroclinicity error which is corrected for below
451         tilde_e3t_a(:,:,:) = MIN( tilde_e3t_a(:,:,:),   rn_zdef_max * e3t_0(:,:,:) )
452         tilde_e3t_a(:,:,:) = MAX( tilde_e3t_a(:,:,:), - rn_zdef_max * e3t_0(:,:,:) )
453
454         !
455         ! "tilda" change in the after scale factor
456         ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
457         DO jk = 1, jpkm1
458            dtilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - tilde_e3t_b(:,:,jk)
459         END DO
460         ! III - Barotropic repartition of the sea surface height over the baroclinic profile
461         ! ==================================================================================
462         ! add ( ssh increment + "baroclinicity error" ) proportionly to e3t(n)
463         ! - ML - baroclinicity error should be better treated in the future
464         !        i.e. locally and not spread over the water column.
465         !        (keep in mind that the idea is to reduce Eulerian velocity as much as possible)
466         zht(:,:) = 0.
467         DO jk = 1, jpkm1
468            zht(:,:)  = zht(:,:) + tilde_e3t_a(:,:,jk) * tmask(:,:,jk)
469         END DO
470         z_scale(:,:) =  - zht(:,:) / ( ht_0(:,:) + sshn(:,:) + 1. - ssmask(:,:) )
471         DO jk = 1, jpkm1
472            dtilde_e3t_a(:,:,jk) = dtilde_e3t_a(:,:,jk) + fse3t_n(:,:,jk) * z_scale(:,:) * tmask(:,:,jk)
473         END DO
474
475      ENDIF
476
477      IF( ln_vvl_ztilde .OR. ln_vvl_layer )  THEN   ! z_tilde or layer coordinate !
478      !                                           ! ---baroclinic part--------- !
479         DO jk = 1, jpkm1
480            fse3t_a(:,:,jk) = fse3t_a(:,:,jk) + dtilde_e3t_a(:,:,jk) * tmask(:,:,jk)
481         END DO
482      ENDIF
483
484      IF( ln_vvl_dbg .AND. .NOT. ll_do_bclinic ) THEN   ! - ML - test: control prints for debuging
485         !
486         IF( lwp ) WRITE(numout, *) 'kt =', kt
487         IF ( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
488            z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( zht(:,:) ) )
489            IF( lk_mpp ) CALL mpp_max( z_tmax )                             ! max over the global domain
490            IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(SUM(tilde_e3t_a))) =', z_tmax
491         END IF
492         !
493         zht(:,:) = 0.0_wp
494         DO jk = 1, jpkm1
495            zht(:,:) = zht(:,:) + fse3t_n(:,:,jk) * tmask(:,:,jk)
496         END DO
497         z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + sshn(:,:) - zht(:,:) ) )
498         IF( lk_mpp ) CALL mpp_max( z_tmax )                                ! max over the global domain
499         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+sshn-SUM(fse3t_n))) =', z_tmax
500         !
501         zht(:,:) = 0.0_wp
502         DO jk = 1, jpkm1
503            zht(:,:) = zht(:,:) + fse3t_a(:,:,jk) * tmask(:,:,jk)
504         END DO
505         z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + ssha(:,:) - zht(:,:) ) )
506         IF( lk_mpp ) CALL mpp_max( z_tmax )                                ! max over the global domain
507         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+ssha-SUM(fse3t_a))) =', z_tmax
508         !
509         zht(:,:) = 0.0_wp
510         DO jk = 1, jpkm1
511            zht(:,:) = zht(:,:) + fse3t_b(:,:,jk) * tmask(:,:,jk)
512         END DO
513         z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + sshb(:,:) - zht(:,:) ) )
514         IF( lk_mpp ) CALL mpp_max( z_tmax )                                ! max over the global domain
515         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+sshb-SUM(fse3t_b))) =', z_tmax
516         !
517         z_tmax = MAXVAL( tmask(:,:,1) *  ABS( sshb(:,:) ) )
518         IF( lk_mpp ) CALL mpp_max( z_tmax )                                ! max over the global domain
519         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(sshb))) =', z_tmax
520         !
521         z_tmax = MAXVAL( tmask(:,:,1) *  ABS( sshn(:,:) ) )
522         IF( lk_mpp ) CALL mpp_max( z_tmax )                                ! max over the global domain
523         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(sshn))) =', z_tmax
524         !
525         z_tmax = MAXVAL( tmask(:,:,1) *  ABS( ssha(:,:) ) )
526         IF( lk_mpp ) CALL mpp_max( z_tmax )                                ! max over the global domain
527         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ssha))) =', z_tmax
528      END IF
529
530      ! *********************************** !
531      ! After scale factors at u- v- points !
532      ! *********************************** !
533
534      CALL dom_vvl_interpol( fse3t_a(:,:,:), fse3u_a(:,:,:), 'U' )
535      CALL dom_vvl_interpol( fse3t_a(:,:,:), fse3v_a(:,:,:), 'V' )
536
537      ! *********************************** !
538      ! After depths at u- v points         !
539      ! *********************************** !
540
541      hu_a(:,:) = 0._wp                        ! Ocean depth at U-points
542      hv_a(:,:) = 0._wp                        ! Ocean depth at V-points
543      DO jk = 1, jpkm1
544         hu_a(:,:) = hu_a(:,:) + fse3u_a(:,:,jk) * umask(:,:,jk)
545         hv_a(:,:) = hv_a(:,:) + fse3v_a(:,:,jk) * vmask(:,:,jk)
546      END DO
547      !                                        ! Inverse of the local depth
548      hur_a(:,:) = 1._wp / ( hu_a(:,:) + 1._wp - umask_i(:,:) ) * umask_i(:,:)
549      hvr_a(:,:) = 1._wp / ( hv_a(:,:) + 1._wp - vmask_i(:,:) ) * vmask_i(:,:)
550
551      CALL wrk_dealloc( jpi, jpj, zht, z_scale, zwu, zwv, zhdiv )
552      CALL wrk_dealloc( jpi, jpj, jpk, ze3t                     )
553
554      IF( nn_timing == 1 )  CALL timing_stop('dom_vvl_sf_nxt')
555
556   END SUBROUTINE dom_vvl_sf_nxt
557
558
559   SUBROUTINE dom_vvl_sf_swp( kt )
560      !!----------------------------------------------------------------------
561      !!                ***  ROUTINE dom_vvl_sf_swp  ***
562      !!                   
563      !! ** Purpose :  compute time filter and swap of scale factors
564      !!               compute all depths and related variables for next time step
565      !!               write outputs and restart file
566      !!
567      !! ** Method  :  - swap of e3t with trick for volume/tracer conservation
568      !!               - reconstruct scale factor at other grid points (interpolate)
569      !!               - recompute depths and water height fields
570      !!
571      !! ** Action  :  - fse3t_(b/n), tilde_e3t_(b/n) and fse3(u/v)_n ready for next time step
572      !!               - Recompute:
573      !!                    fse3(u/v)_b       
574      !!                    fse3w_n           
575      !!                    fse3(u/v)w_b     
576      !!                    fse3(u/v)w_n     
577      !!                    fsdept_n, fsdepw_n  and fsde3w_n
578      !!                    h(u/v) and h(u/v)r
579      !!
580      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
581      !!              Leclair, M., and G. Madec, 2011, Ocean Modelling.
582      !!----------------------------------------------------------------------
583      !! * Arguments
584      INTEGER, INTENT( in )               :: kt       ! time step
585      !! * Local declarations
586      INTEGER                             :: ji,jj,jk       ! dummy loop indices
587      REAL(wp)                            :: zcoef
588      !!----------------------------------------------------------------------
589
590      IF( nn_timing == 1 )  CALL timing_start('dom_vvl_sf_swp')
591      !
592      IF( kt == nit000 )   THEN
593         IF(lwp) WRITE(numout,*)
594         IF(lwp) WRITE(numout,*) 'dom_vvl_sf_swp : - time filter and swap of scale factors'
595         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~   - interpolate scale factors and compute depths for next time step'
596      ENDIF
597
598      !
599      ! Time filter and swap of scale factors
600      ! =====================================
601      ! - ML - fse3(t/u/v)_b are allready computed in dynnxt.
602      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
603         IF( neuler == 0 .AND. kt == nit000 ) THEN
604            tilde_e3t_b(:,:,:) = tilde_e3t_n(:,:,:)
605         ELSE
606            tilde_e3t_b(:,:,:) = tilde_e3t_n(:,:,:) & 
607            &         + atfp * ( tilde_e3t_b(:,:,:) - 2.0_wp * tilde_e3t_n(:,:,:) + tilde_e3t_a(:,:,:) )
608         ENDIF
609         tilde_e3t_n(:,:,:) = tilde_e3t_a(:,:,:)
610      ENDIF
611      fsdept_b(:,:,:) = fsdept_n(:,:,:)
612      fsdepw_b(:,:,:) = fsdepw_n(:,:,:)
613
614      fse3t_n(:,:,:) = fse3t_a(:,:,:)
615      fse3u_n(:,:,:) = fse3u_a(:,:,:)
616      fse3v_n(:,:,:) = fse3v_a(:,:,:)
617
618      ! Compute all missing vertical scale factor and depths
619      ! ====================================================
620      ! Horizontal scale factor interpolations
621      ! --------------------------------------
622      ! - ML - fse3u_b and fse3v_b are allready computed in dynnxt
623      ! - JC - hu_b, hv_b, hur_b, hvr_b also
624      CALL dom_vvl_interpol( fse3u_n(:,:,:), fse3f_n (:,:,:), 'F'  )
625      ! Vertical scale factor interpolations
626      ! ------------------------------------
627      CALL dom_vvl_interpol( fse3t_n(:,:,:), fse3w_n (:,:,:), 'W'  )
628      CALL dom_vvl_interpol( fse3u_n(:,:,:), fse3uw_n(:,:,:), 'UW' )
629      CALL dom_vvl_interpol( fse3v_n(:,:,:), fse3vw_n(:,:,:), 'VW' )
630      CALL dom_vvl_interpol( fse3t_b(:,:,:), fse3w_b (:,:,:), 'W'  )
631      CALL dom_vvl_interpol( fse3u_b(:,:,:), fse3uw_b(:,:,:), 'UW' )
632      CALL dom_vvl_interpol( fse3v_b(:,:,:), fse3vw_b(:,:,:), 'VW' )
633      ! t- and w- points depth
634      ! ----------------------
635      ! set the isf depth as it is in the initial step
636      fsdept_n(:,:,1) = 0.5_wp * fse3w_n(:,:,1)
637      fsdepw_n(:,:,1) = 0.0_wp
638      fsde3w_n(:,:,1) = fsdept_n(:,:,1) - sshn(:,:)
639
640      DO jk = 2, jpk
641         DO jj = 1,jpj
642            DO ji = 1,jpi
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               fsdepw_n(ji,jj,jk) = fsdepw_n(ji,jj,jk-1) + fse3t_n(ji,jj,jk-1)
647               fsdept_n(ji,jj,jk) =      zcoef  * ( fsdepw_n(ji,jj,jk  ) + 0.5 * fse3w_n(ji,jj,jk))  &
648                   &                + (1-zcoef) * ( fsdept_n(ji,jj,jk-1) +       fse3w_n(ji,jj,jk)) 
649               fsde3w_n(ji,jj,jk) = fsdept_n(ji,jj,jk) - sshn(ji,jj)
650            END DO
651         END DO
652      END DO
653
654      ! Local depth and Inverse of the local depth of the water column at u- and v- points
655      ! ----------------------------------------------------------------------------------
656      hu (:,:) = hu_a (:,:)
657      hv (:,:) = hv_a (:,:)
658
659      ! Inverse of the local depth
660      hur(:,:) = hur_a(:,:)
661      hvr(:,:) = hvr_a(:,:)
662
663      ! Local depth of the water column at t- points
664      ! --------------------------------------------
665      ht(:,:) = 0.
666      DO jk = 1, jpkm1
667         ht(:,:) = ht(:,:) + fse3t_n(:,:,jk) * tmask(:,:,jk)
668      END DO
669
670      ! write restart file
671      ! ==================
672      IF( lrst_oce ) CALL dom_vvl_rst( kt, 'WRITE' )
673      !
674      IF( nn_timing == 1 )  CALL timing_stop('dom_vvl_sf_swp')
675
676   END SUBROUTINE dom_vvl_sf_swp
677
678
679   SUBROUTINE dom_vvl_interpol( pe3_in, pe3_out, pout )
680      !!---------------------------------------------------------------------
681      !!                  ***  ROUTINE dom_vvl__interpol  ***
682      !!
683      !! ** Purpose :   interpolate scale factors from one grid point to another
684      !!
685      !! ** Method  :   e3_out = e3_0 + interpolation(e3_in - e3_0)
686      !!                - horizontal interpolation: grid cell surface averaging
687      !!                - vertical interpolation: simple averaging
688      !!----------------------------------------------------------------------
689      !! * Arguments
690      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( in    ) ::  pe3_in     ! input e3 to be interpolated
691      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( inout ) ::  pe3_out    ! output interpolated e3
692      CHARACTER(LEN=*), INTENT( in )                    ::  pout       ! grid point of out scale factors
693      !                                                                !   =  'U', 'V', 'W, 'F', 'UW' or 'VW'
694      !! * Local declarations
695      INTEGER ::   ji, jj, jk                                          ! dummy loop indices
696      LOGICAL ::   l_is_orca                                           ! local logical
697      !!----------------------------------------------------------------------
698      IF( nn_timing == 1 )  CALL timing_start('dom_vvl_interpol')
699         !
700      l_is_orca = .FALSE.
701      IF( cp_cfg == "orca" .AND. jp_cfg == 2 ) l_is_orca = .TRUE.      ! ORCA R2 configuration - will need to correct some locations
702
703      SELECT CASE ( pout )
704         !               ! ------------------------------------- !
705      CASE( 'U' )        ! interpolation from T-point to U-point !
706         !               ! ------------------------------------- !
707         ! horizontal surface weighted interpolation
708         DO jk = 1, jpk
709            DO jj = 1, jpjm1
710               DO ji = 1, fs_jpim1   ! vector opt.
711                  pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk) * r1_e12u(ji,jj)                                   &
712                     &                       * (   e12t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) )     &
713                     &                           + e12t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) )
714               END DO
715            END DO
716         END DO
717         !
718         IF( l_is_orca ) CALL dom_vvl_orca_fix( pe3_in, pe3_out, pout )
719         ! boundary conditions
720         CALL lbc_lnk( pe3_out(:,:,:), 'U', 1._wp )
721         pe3_out(:,:,:) = pe3_out(:,:,:) + e3u_0(:,:,:)
722         !               ! ------------------------------------- !
723      CASE( 'V' )        ! interpolation from T-point to V-point !
724         !               ! ------------------------------------- !
725         ! horizontal surface weighted interpolation
726         DO jk = 1, jpk
727            DO jj = 1, jpjm1
728               DO ji = 1, fs_jpim1   ! vector opt.
729                  pe3_out(ji,jj,jk) = 0.5_wp * vmask(ji,jj,jk) * r1_e12v(ji,jj)                                   &
730                     &                       * (   e12t(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3t_0(ji,jj  ,jk) )     &
731                     &                           + e12t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) )
732               END DO
733            END DO
734         END DO
735         !
736         IF( l_is_orca ) CALL dom_vvl_orca_fix( pe3_in, pe3_out, pout )
737         ! boundary conditions
738         CALL lbc_lnk( pe3_out(:,:,:), 'V', 1._wp )
739         pe3_out(:,:,:) = pe3_out(:,:,:) + e3v_0(:,:,:)
740         !               ! ------------------------------------- !
741      CASE( 'F' )        ! interpolation from U-point to F-point !
742         !               ! ------------------------------------- !
743         ! horizontal surface weighted interpolation
744         DO jk = 1, jpk
745            DO jj = 1, jpjm1
746               DO ji = 1, fs_jpim1   ! vector opt.
747                  pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk) * umask(ji,jj+1,jk) * r1_e12f(ji,jj)               &
748                     &                       * (   e12u(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3u_0(ji,jj  ,jk) )     &
749                     &                           + e12u(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3u_0(ji,jj+1,jk) ) )
750               END DO
751            END DO
752         END DO
753         !
754         IF( l_is_orca ) CALL dom_vvl_orca_fix( pe3_in, pe3_out, pout )
755         ! boundary conditions
756         CALL lbc_lnk( pe3_out(:,:,:), 'F', 1._wp )
757         pe3_out(:,:,:) = pe3_out(:,:,:) + e3f_0(:,:,:)
758         !               ! ------------------------------------- !
759      CASE( 'W' )        ! interpolation from T-point to W-point !
760         !               ! ------------------------------------- !
761         ! vertical simple interpolation
762         pe3_out(:,:,1) = e3w_0(:,:,1) + pe3_in(:,:,1) - e3t_0(:,:,1)
763         ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without indirect adressing
764         DO jk = 2, jpk
765            pe3_out(:,:,jk) = e3w_0(:,:,jk) + ( 1.0_wp - 0.5_wp * tmask(:,:,jk) ) * ( pe3_in(:,:,jk-1) - e3t_0(:,:,jk-1) )   &
766               &                            +            0.5_wp * tmask(:,:,jk)   * ( pe3_in(:,:,jk  ) - e3t_0(:,:,jk  ) )
767         END DO
768         !               ! -------------------------------------- !
769      CASE( 'UW' )       ! interpolation from U-point to UW-point !
770         !               ! -------------------------------------- !
771         ! vertical simple interpolation
772         pe3_out(:,:,1) = e3uw_0(:,:,1) + pe3_in(:,:,1) - e3u_0(:,:,1)
773         ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without indirect adressing
774         DO jk = 2, jpk
775            pe3_out(:,:,jk) = e3uw_0(:,:,jk) + ( 1.0_wp - 0.5_wp * umask(:,:,jk) ) * ( pe3_in(:,:,jk-1) - e3u_0(:,:,jk-1) )   &
776               &                             +            0.5_wp * umask(:,:,jk)   * ( pe3_in(:,:,jk  ) - e3u_0(:,:,jk  ) )
777         END DO
778         !               ! -------------------------------------- !
779      CASE( 'VW' )       ! interpolation from V-point to VW-point !
780         !               ! -------------------------------------- !
781         ! vertical simple interpolation
782         pe3_out(:,:,1) = e3vw_0(:,:,1) + pe3_in(:,:,1) - e3v_0(:,:,1)
783         ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without indirect adressing
784         DO jk = 2, jpk
785            pe3_out(:,:,jk) = e3vw_0(:,:,jk) + ( 1.0_wp - 0.5_wp * vmask(:,:,jk) ) * ( pe3_in(:,:,jk-1) - e3v_0(:,:,jk-1) )   &
786               &                             +            0.5_wp * vmask(:,:,jk)   * ( pe3_in(:,:,jk  ) - e3v_0(:,:,jk  ) )
787         END DO
788      END SELECT
789      !
790
791      IF( nn_timing == 1 )  CALL timing_stop('dom_vvl_interpol')
792
793   END SUBROUTINE dom_vvl_interpol
794
795   SUBROUTINE dom_vvl_rst( kt, cdrw )
796      !!---------------------------------------------------------------------
797      !!                   ***  ROUTINE dom_vvl_rst  ***
798      !!                     
799      !! ** Purpose :   Read or write VVL file in restart file
800      !!
801      !! ** Method  :   use of IOM library
802      !!                if the restart does not contain vertical scale factors,
803      !!                they are set to the _0 values
804      !!                if the restart does not contain vertical scale factors increments (z_tilde),
805      !!                they are set to 0.
806      !!----------------------------------------------------------------------
807      !! * Arguments
808      INTEGER         , INTENT(in) ::   kt     ! ocean time-step
809      CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
810      !! * Local declarations
811      INTEGER ::   jk
812      INTEGER ::   id1, id2, id3, id4, id5     ! local integers
813      !!----------------------------------------------------------------------
814      !
815      IF( nn_timing == 1 )  CALL timing_start('dom_vvl_rst')
816      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
817         !                                   ! ===============
818         IF( ln_rstart ) THEN                   !* Read the restart file
819            CALL rst_read_open                  !  open the restart file if necessary
820            CALL iom_get( numror, jpdom_autoglo, 'sshn'   , sshn, lrxios = lxios_read    )
821            !
822            id1 = iom_varid( numror, 'fse3t_b', ldstop = .FALSE. )
823            id2 = iom_varid( numror, 'fse3t_n', ldstop = .FALSE. )
824            id3 = iom_varid( numror, 'tilde_e3t_b', ldstop = .FALSE. )
825            id4 = iom_varid( numror, 'tilde_e3t_n', ldstop = .FALSE. )
826            id5 = iom_varid( numror, 'hdiv_lf', ldstop = .FALSE. )
827            !                             ! --------- !
828            !                             ! all cases !
829            !                             ! --------- !
830            IF( MIN( id1, id2 ) > 0 ) THEN       ! all required arrays exist
831               CALL iom_get( numror, jpdom_autoglo, 'fse3t_b', fse3t_b(:,:,:), lrxios = lxios_read )
832               CALL iom_get( numror, jpdom_autoglo, 'fse3t_n', fse3t_n(:,:,:), lrxios = lxios_read )
833               ! needed to restart if land processor not computed
834               IF(lwp) write(numout,*) 'dom_vvl_rst : fse3t_b and fse3t_n found in restart files'
835               WHERE ( tmask(:,:,:) == 0.0_wp ) 
836                  fse3t_n(:,:,:) = e3t_0(:,:,:)
837                  fse3t_b(:,:,:) = e3t_0(:,:,:)
838               END WHERE
839               IF( neuler == 0 ) THEN
840                  fse3t_b(:,:,:) = fse3t_n(:,:,:)
841               ENDIF
842            ELSE IF( id1 > 0 ) THEN
843               IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : fse3t_n not found in restart files'
844               IF(lwp) write(numout,*) 'fse3t_n set equal to fse3t_b.'
845               IF(lwp) write(numout,*) 'neuler is forced to 0'
846               CALL iom_get( numror, jpdom_autoglo, 'fse3t_b', fse3t_b(:,:,:), lrxios = lxios_read )
847               fse3t_n(:,:,:) = fse3t_b(:,:,:)
848               neuler = 0
849            ELSE IF( id2 > 0 ) THEN
850               IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : fse3t_b not found in restart files'
851               IF(lwp) write(numout,*) 'fse3t_b set equal to fse3t_n.'
852               IF(lwp) write(numout,*) 'neuler is forced to 0'
853               CALL iom_get( numror, jpdom_autoglo, 'fse3t_n', fse3t_n(:,:,:), lrxios = lxios_read )
854               fse3t_b(:,:,:) = fse3t_n(:,:,:)
855               neuler = 0
856            ELSE
857               IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : fse3t_n not found in restart file'
858               IF(lwp) write(numout,*) 'Compute scale factor from sshn'
859               IF(lwp) write(numout,*) 'neuler is forced to 0'
860               DO jk=1,jpk
861                  fse3t_n(:,:,jk) =  e3t_0(:,:,jk) * ( ht_0(:,:) + sshn(:,:) ) &
862                      &                            / ( ht_0(:,:) + 1._wp - ssmask(:,:) ) * tmask(:,:,jk) &
863                      &            + e3t_0(:,:,jk) * (1._wp -tmask(:,:,jk))
864               END DO
865               fse3t_b(:,:,:) = fse3t_n(:,:,:)
866               neuler = 0
867            ENDIF
868            !                             ! ----------- !
869            IF( ln_vvl_zstar ) THEN       ! z_star case !
870               !                          ! ----------- !
871               IF( MIN( id3, id4 ) > 0 ) THEN
872                  CALL ctl_stop( 'dom_vvl_rst: z_star cannot restart from a z_tilde or layer run' )
873               ENDIF
874               !                          ! ----------------------- !
875            ELSE                          ! z_tilde and layer cases !
876               !                          ! ----------------------- !
877               IF( MIN( id3, id4 ) > 0 ) THEN  ! all required arrays exist
878                  CALL iom_get( numror, jpdom_autoglo, 'tilde_e3t_b', tilde_e3t_b(:,:,:), lrxios = lxios_read )
879                  CALL iom_get( numror, jpdom_autoglo, 'tilde_e3t_n', tilde_e3t_n(:,:,:), lrxios = lxios_read )
880               ELSE                            ! one at least array is missing
881                  tilde_e3t_b(:,:,:) = 0.0_wp
882                  tilde_e3t_n(:,:,:) = 0.0_wp
883               ENDIF
884               !                          ! ------------ !
885               IF( ln_vvl_ztilde ) THEN   ! z_tilde case !
886                  !                       ! ------------ !
887                  IF( id5 > 0 ) THEN  ! required array exists
888                     CALL iom_get( numror, jpdom_autoglo, 'hdiv_lf', hdiv_lf(:,:,:), lrxios = lxios_read )
889                  ELSE                ! array is missing
890                     hdiv_lf(:,:,:) = 0.0_wp
891                  ENDIF
892               ENDIF
893            ENDIF
894            !
895         ELSE                                   !* Initialize at "rest"
896            fse3t_b(:,:,:) = e3t_0(:,:,:)
897            fse3t_n(:,:,:) = e3t_0(:,:,:)
898            sshn(:,:) = 0.0_wp
899            IF( ln_vvl_ztilde .OR. ln_vvl_layer) THEN
900               tilde_e3t_b(:,:,:) = 0.0_wp
901               tilde_e3t_n(:,:,:) = 0.0_wp
902               IF( ln_vvl_ztilde ) hdiv_lf(:,:,:) = 0.0_wp
903            END IF
904         ENDIF
905
906      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
907         !                                   ! ===================
908         IF(lwp) WRITE(numout,*) '---- dom_vvl_rst ----'
909         !                                           ! --------- !
910         !                                           ! all cases !
911         !                                           ! --------- !
912         CALL iom_rstput( kt, nitrst, numrow, 'fse3t_b', fse3t_b(:,:,:) )
913         CALL iom_rstput( kt, nitrst, numrow, 'fse3t_n', fse3t_n(:,:,:) )
914         !                                           ! ----------------------- !
915         IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN  ! z_tilde and layer cases !
916            !                                        ! ----------------------- !
917            CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_b', tilde_e3t_b(:,:,:) )
918            CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_n', tilde_e3t_n(:,:,:) )
919         END IF
920         !                                           ! -------------!   
921         IF( ln_vvl_ztilde ) THEN                    ! z_tilde case !
922            !                                        ! ------------ !
923            CALL iom_rstput( kt, nitrst, numrow, 'hdiv_lf', hdiv_lf(:,:,:) )
924         ENDIF
925
926      ENDIF
927      IF( nn_timing == 1 )  CALL timing_stop('dom_vvl_rst')
928
929   END SUBROUTINE dom_vvl_rst
930
931
932   SUBROUTINE dom_vvl_ctl
933      !!---------------------------------------------------------------------
934      !!                  ***  ROUTINE dom_vvl_ctl  ***
935      !!               
936      !! ** Purpose :   Control the consistency between namelist options
937      !!                for vertical coordinate
938      !!----------------------------------------------------------------------
939      INTEGER ::   ioptio
940      INTEGER ::   ios
941
942      NAMELIST/nam_vvl/ ln_vvl_zstar, ln_vvl_ztilde, ln_vvl_layer, ln_vvl_ztilde_as_zstar, &
943                      & ln_vvl_zstar_at_eqtor      , rn_ahe3     , rn_rst_e3t            , &
944                      & rn_lf_cutoff               , rn_zdef_max , ln_vvl_dbg                ! not yet implemented: ln_vvl_kepe
945      !!----------------------------------------------------------------------
946
947      REWIND( numnam_ref )              ! Namelist nam_vvl in reference namelist :
948      READ  ( numnam_ref, nam_vvl, IOSTAT = ios, ERR = 901)
949901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nam_vvl in reference namelist', lwp )
950
951      REWIND( numnam_cfg )              ! Namelist nam_vvl in configuration namelist : Parameters of the run
952      READ  ( numnam_cfg, nam_vvl, IOSTAT = ios, ERR = 902 )
953902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nam_vvl in configuration namelist', lwp )
954      IF(lwm) WRITE ( numond, nam_vvl )
955
956      IF(lwp) THEN                    ! Namelist print
957         WRITE(numout,*)
958         WRITE(numout,*) 'dom_vvl_ctl : choice/control of the variable vertical coordinate'
959         WRITE(numout,*) '~~~~~~~~~~~'
960         WRITE(numout,*) '           Namelist nam_vvl : chose a vertical coordinate'
961         WRITE(numout,*) '              zstar                      ln_vvl_zstar   = ', ln_vvl_zstar
962         WRITE(numout,*) '              ztilde                     ln_vvl_ztilde  = ', ln_vvl_ztilde
963         WRITE(numout,*) '              layer                      ln_vvl_layer   = ', ln_vvl_layer
964         WRITE(numout,*) '              ztilde as zstar   ln_vvl_ztilde_as_zstar  = ', ln_vvl_ztilde_as_zstar
965         WRITE(numout,*) '      ztilde near the equator    ln_vvl_zstar_at_eqtor  = ', ln_vvl_zstar_at_eqtor
966         ! WRITE(numout,*) '           Namelist nam_vvl : chose kinetic-to-potential energy conservation'
967         ! WRITE(numout,*) '                                         ln_vvl_kepe    = ', ln_vvl_kepe
968         WRITE(numout,*) '           Namelist nam_vvl : thickness diffusion coefficient'
969         WRITE(numout,*) '                                         rn_ahe3        = ', rn_ahe3
970         WRITE(numout,*) '           Namelist nam_vvl : maximum e3t deformation fractional change'
971         WRITE(numout,*) '                                         rn_zdef_max    = ', rn_zdef_max
972         IF( ln_vvl_ztilde_as_zstar ) THEN
973            WRITE(numout,*) '           ztilde running in zstar emulation mode; '
974            WRITE(numout,*) '           ignoring namelist timescale parameters and using:'
975            WRITE(numout,*) '                 hard-wired : z-tilde to zstar restoration timescale (days)'
976            WRITE(numout,*) '                                         rn_rst_e3t     =    0.0'
977            WRITE(numout,*) '                 hard-wired : z-tilde cutoff frequency of low-pass filter (days)'
978            WRITE(numout,*) '                                         rn_lf_cutoff   =    1.0/rdt'
979         ELSE
980            WRITE(numout,*) '           Namelist nam_vvl : z-tilde to zstar restoration timescale (days)'
981            WRITE(numout,*) '                                         rn_rst_e3t     = ', rn_rst_e3t
982            WRITE(numout,*) '           Namelist nam_vvl : z-tilde cutoff frequency of low-pass filter (days)'
983            WRITE(numout,*) '                                         rn_lf_cutoff   = ', rn_lf_cutoff
984         ENDIF
985         WRITE(numout,*) '           Namelist nam_vvl : debug prints'
986         WRITE(numout,*) '                                         ln_vvl_dbg     = ', ln_vvl_dbg
987      ENDIF
988
989      ioptio = 0                      ! Parameter control
990      IF( ln_vvl_ztilde_as_zstar ) ln_vvl_ztilde = .true.
991      IF( ln_vvl_zstar           )        ioptio = ioptio + 1
992      IF( ln_vvl_ztilde          )        ioptio = ioptio + 1
993      IF( ln_vvl_layer           )        ioptio = ioptio + 1
994
995      IF( ioptio /= 1 )   CALL ctl_stop( 'Choose ONE vertical coordinate in namelist nam_vvl' )
996      IF( .NOT. ln_vvl_zstar .AND. nn_isf .NE. 0) CALL ctl_stop( 'Only vvl_zstar has been tested with ice shelf cavity' )
997
998      IF(lwp) THEN                   ! Print the choice
999         WRITE(numout,*)
1000         IF( ln_vvl_zstar           ) WRITE(numout,*) '              zstar vertical coordinate is used'
1001         IF( ln_vvl_ztilde          ) WRITE(numout,*) '              ztilde vertical coordinate is used'
1002         IF( ln_vvl_layer           ) WRITE(numout,*) '              layer vertical coordinate is used'
1003         IF( ln_vvl_ztilde_as_zstar ) WRITE(numout,*) '              to emulate a zstar coordinate'
1004         ! - ML - Option not developed yet
1005         ! IF(       ln_vvl_kepe ) WRITE(numout,*) '              kinetic to potential energy transfer : option used'
1006         ! IF( .NOT. ln_vvl_kepe ) WRITE(numout,*) '              kinetic to potential energy transfer : option not used'
1007      ENDIF
1008
1009#if defined key_agrif
1010      IF (.NOT.Agrif_Root()) CALL ctl_stop( 'AGRIF not implemented with non-linear free surface (key_vvl)' )
1011#endif
1012
1013   END SUBROUTINE dom_vvl_ctl
1014
1015   SUBROUTINE dom_vvl_orca_fix( pe3_in, pe3_out, pout )
1016      !!---------------------------------------------------------------------
1017      !!                   ***  ROUTINE dom_vvl_orca_fix  ***
1018      !!                     
1019      !! ** Purpose :   Correct surface weighted, horizontally interpolated,
1020      !!                scale factors at locations that have been individually
1021      !!                modified in domhgr. Such modifications break the
1022      !!                relationship between e12t and e1u*e2u etc.
1023      !!                Recompute some scale factors ignoring the modified metric.
1024      !!----------------------------------------------------------------------
1025      !! * Arguments
1026      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( in    ) ::  pe3_in     ! input e3 to be interpolated
1027      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT( inout ) ::  pe3_out    ! output interpolated e3
1028      CHARACTER(LEN=*), INTENT( in )                    ::  pout       ! grid point of out scale factors
1029      !                                                                !   =  'U', 'V', 'W, 'F', 'UW' or 'VW'
1030      !! * Local declarations
1031      INTEGER ::   ji, jj, jk                                          ! dummy loop indices
1032      INTEGER ::   ij0, ij1, ii0, ii1                                  ! dummy loop indices
1033      INTEGER ::   isrow                                               ! index for ORCA1 starting row
1034      !! acc
1035      !! Hmm with the time splitting these "fixes" seem to do more harm than good. Temporarily disabled for
1036      !! the ORCA2 tests (by changing jp_cfg test from 2 to 3) pending further investigations
1037      !!
1038      !                                                ! =====================
1039      IF( cp_cfg == "orca" .AND. jp_cfg == 3 ) THEN    ! ORCA R2 configuration
1040         !                                             ! =====================
1041      !! acc
1042         IF( nn_cla == 0 ) THEN
1043            !
1044            ii0 = 139   ;   ii1 = 140        ! Gibraltar Strait (e2u was modified)
1045            ij0 = 102   ;   ij1 = 102
1046            DO jk = 1, jpkm1
1047               DO jj = mj0(ij0), mj1(ij1)
1048                  DO ji = mi0(ii0), mi1(ii1)
1049                     SELECT CASE ( pout )
1050                     CASE( 'U' )
1051                        pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk)                                        &
1052                       &                    * (   e1t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) ) &
1053                       &                    +     e1t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) &
1054                       &                      ) / e1u(ji,jj)   +   e3u_0(ji,jj,jk)
1055                     CASE( 'F' )
1056                        pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk) * umask(ji,jj+1,jk)                    &
1057                       &                    * (   e1u(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3u_0(ji  ,jj,jk) ) &
1058                       &                    +     e1u(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3u_0(ji+1,jj,jk) ) &
1059                       &                      ) / e1f(ji,jj)   +   e3f_0(ji,jj,jk)
1060                     END SELECT
1061                  END DO
1062               END DO
1063            END DO
1064            !
1065            ii0 = 160   ;   ii1 = 160        ! Bab el Mandeb (e2u and e1v were modified)
1066            ij0 =  88   ;   ij1 =  88
1067            DO jk = 1, jpkm1
1068               DO jj = mj0(ij0), mj1(ij1)
1069                  DO ji = mi0(ii0), mi1(ii1)
1070                     SELECT CASE ( pout )
1071                     CASE( 'U' )
1072                        pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk)                                        &
1073                       &                    * (   e1t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) ) &
1074                       &                    +     e1t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) &
1075                       &                      ) / e1u(ji,jj)   +   e3u_0(ji,jj,jk)
1076                     CASE( 'V' )
1077                        pe3_out(ji,jj,jk) = 0.5_wp * vmask(ji,jj,jk)                                        &
1078                       &                    * (   e2t(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3t_0(ji,jj  ,jk) ) &
1079                       &                    +     e2t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) &
1080                       &                      ) / e2v(ji,jj)   +   e3v_0(ji,jj,jk)
1081                     CASE( 'F' )
1082                        pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk) * umask(ji,jj+1,jk)                    &
1083                       &                    * (   e1u(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3u_0(ji  ,jj,jk) ) &
1084                       &                    +     e1u(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3u_0(ji+1,jj,jk) ) &
1085                       &                      ) / e1f(ji,jj)   +   e3f_0(ji,jj,jk)
1086                     END SELECT
1087                  END DO
1088               END DO
1089            END DO
1090         ENDIF
1091
1092         ii0 = 145   ;   ii1 = 146        ! Danish Straits (e2u was modified)
1093         ij0 = 116   ;   ij1 = 116
1094         DO jk = 1, jpkm1
1095            DO jj = mj0(ij0), mj1(ij1)
1096               DO ji = mi0(ii0), mi1(ii1)
1097                  SELECT CASE ( pout )
1098                  CASE( 'U' )
1099                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk)                                        &
1100                    &                    * (   e1t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) ) &
1101                    &                    +     e1t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) &
1102                    &                      ) / e1u(ji,jj)   +   e3u_0(ji,jj,jk)
1103                  CASE( 'F' )
1104                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk) * umask(ji,jj+1,jk)                    &
1105                    &                    * (   e1u(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3u_0(ji  ,jj,jk) ) &
1106                    &                    +     e1u(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3u_0(ji+1,jj,jk) ) &
1107                    &                      ) / e1f(ji,jj)   +   e3f_0(ji,jj,jk)
1108                  END SELECT
1109               END DO
1110            END DO
1111         END DO
1112      ENDIF
1113      !
1114         !                                             ! =====================
1115      IF( cp_cfg == "orca" .AND. jp_cfg == 1 ) THEN    ! ORCA R1 configuration
1116         !                                             ! =====================
1117         ! This dirty section will be suppressed by simplification process:
1118         ! all this will come back in input files
1119         ! Currently these hard-wired indices relate to configuration with
1120         ! extend grid (jpjglo=332)
1121         ! which had a grid-size of 362x292.
1122         isrow = 332 - jpjglo
1123         !
1124         ii0 = 282           ;   ii1 = 283        ! Gibraltar Strait (e2u was modified)
1125         ij0 = 241 - isrow   ;   ij1 = 241 - isrow
1126         DO jk = 1, jpkm1
1127            DO jj = mj0(ij0), mj1(ij1)
1128               DO ji = mi0(ii0), mi1(ii1)
1129                  SELECT CASE ( pout )
1130                  CASE( 'U' )
1131                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk)                                        &
1132                    &                    * (   e1t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) ) &
1133                    &                    +     e1t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) &
1134                    &                      ) / e1u(ji,jj)   +   e3u_0(ji,jj,jk)
1135                  CASE( 'F' )
1136                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk) * umask(ji,jj+1,jk)                    &
1137                    &                    * (   e1u(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3u_0(ji  ,jj,jk) ) &
1138                    &                    +     e1u(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3u_0(ji+1,jj,jk) ) &
1139                    &                      ) / e1f(ji,jj)   +   e3f_0(ji,jj,jk)
1140                  END SELECT
1141               END DO
1142            END DO
1143         END DO
1144         !
1145         ii0 = 314           ;   ii1 = 315        ! Bhosporus Strait (e2u was modified)
1146         ij0 = 248 - isrow   ;   ij1 = 248 - isrow
1147         DO jk = 1, jpkm1
1148            DO jj = mj0(ij0), mj1(ij1)
1149               DO ji = mi0(ii0), mi1(ii1)
1150                  SELECT CASE ( pout )
1151                  CASE( 'U' )
1152                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk)                                        & 
1153                    &                    * (   e1t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) ) &
1154                    &                    +     e1t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) &
1155                    &                      ) / e1u(ji,jj)   +   e3u_0(ji,jj,jk)
1156                  CASE( 'F' )
1157                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk) * umask(ji,jj+1,jk)                    & 
1158                    &                    * (   e1u(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3u_0(ji  ,jj,jk) ) &
1159                    &                    +     e1u(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3u_0(ji+1,jj,jk) ) &
1160                    &                      ) / e1f(ji,jj)   +   e3f_0(ji,jj,jk)
1161                  END SELECT
1162               END DO
1163            END DO
1164         END DO
1165         !
1166         ii0 =  44           ;   ii1 =  44        ! Lombok Strait (e1v was modified)
1167         ij0 = 164 - isrow   ;   ij1 = 165 - isrow
1168         DO jk = 1, jpkm1
1169            DO jj = mj0(ij0), mj1(ij1)
1170               DO ji = mi0(ii0), mi1(ii1)
1171                  SELECT CASE ( pout )
1172                  CASE( 'V' )
1173                     pe3_out(ji,jj,jk) = 0.5_wp * vmask(ji,jj,jk)                                        &
1174                    &                    * (   e2t(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3t_0(ji,jj  ,jk) ) &
1175                    &                    +     e2t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) &
1176                    &                      ) / e2v(ji,jj)   +   e3v_0(ji,jj,jk)
1177                  END SELECT
1178               END DO
1179            END DO
1180         END DO
1181         !
1182         ii0 =  48           ;   ii1 =  48        ! Sumba Strait (e1v was modified) [closed from bathy_11 on]
1183         ij0 = 164 - isrow   ;   ij1 = 165 - isrow
1184         DO jk = 1, jpkm1
1185            DO jj = mj0(ij0), mj1(ij1)
1186               DO ji = mi0(ii0), mi1(ii1)
1187                  SELECT CASE ( pout )
1188                  CASE( 'V' )
1189                     pe3_out(ji,jj,jk) = 0.5_wp * vmask(ji,jj,jk)                                        &
1190                    &                    * (   e2t(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3t_0(ji,jj  ,jk) ) &
1191                    &                    +     e2t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) &
1192                    &                      ) / e2v(ji,jj)   +   e3v_0(ji,jj,jk)
1193                  END SELECT
1194               END DO
1195            END DO
1196         END DO
1197         !
1198         ii0 =  53          ;   ii1 =  53        ! Ombai Strait (e1v was modified)
1199         ij0 = 164 - isrow  ;   ij1 = 165  - isrow 
1200         DO jk = 1, jpkm1
1201            DO jj = mj0(ij0), mj1(ij1)
1202               DO ji = mi0(ii0), mi1(ii1)
1203                  SELECT CASE ( pout )
1204                  CASE( 'V' )
1205                     pe3_out(ji,jj,jk) = 0.5_wp * vmask(ji,jj,jk)                                        &
1206                    &                    * (   e2t(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3t_0(ji,jj  ,jk) ) &
1207                    &                    +     e2t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) &
1208                    &                      ) / e2v(ji,jj)   +   e3v_0(ji,jj,jk)
1209                  END SELECT
1210               END DO
1211            END DO
1212         END DO
1213         !
1214         ii0 =  56            ;   ii1 =  56        ! Timor Passage (e1v was modified)
1215         ij0 = 164 - isrow    ;   ij1 = 165  - isrow 
1216         DO jk = 1, jpkm1
1217            DO jj = mj0(ij0), mj1(ij1)
1218               DO ji = mi0(ii0), mi1(ii1)
1219                  SELECT CASE ( pout )
1220                  CASE( 'V' )
1221                     pe3_out(ji,jj,jk) = 0.5_wp * vmask(ji,jj,jk)                                        &
1222                    &                    * (   e2t(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3t_0(ji,jj  ,jk) ) &
1223                    &                    +     e2t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) &
1224                    &                      ) / e2v(ji,jj)   +   e3v_0(ji,jj,jk)
1225                  END SELECT
1226               END DO
1227            END DO
1228         END DO
1229         !
1230         ii0 =  55            ;   ii1 =  55        ! West Halmahera Strait (e1v was modified)
1231         ij0 = 181 - isrow    ;   ij1 = 182 - isrow 
1232         DO jk = 1, jpkm1
1233            DO jj = mj0(ij0), mj1(ij1)
1234               DO ji = mi0(ii0), mi1(ii1)
1235                  SELECT CASE ( pout )
1236                  CASE( 'V' )
1237                     pe3_out(ji,jj,jk) = 0.5_wp * vmask(ji,jj,jk)                                        &
1238                    &                    * (   e2t(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3t_0(ji,jj  ,jk) ) &
1239                    &                    +     e2t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) &
1240                    &                      ) / e2v(ji,jj)   +   e3v_0(ji,jj,jk)
1241                  END SELECT
1242               END DO
1243            END DO
1244         END DO
1245         !
1246         ii0 =  58            ;   ii1 =  58        ! East Halmahera Strait (e1v was modified)
1247         ij0 = 181 - isrow    ;   ij1 = 182 - isrow 
1248         DO jk = 1, jpkm1
1249            DO jj = mj0(ij0), mj1(ij1)
1250               DO ji = mi0(ii0), mi1(ii1)
1251                  SELECT CASE ( pout )
1252                  CASE( 'V' )
1253                     pe3_out(ji,jj,jk) = 0.5_wp * vmask(ji,jj,jk)                                        &
1254                    &                    * (   e2t(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3t_0(ji,jj  ,jk) ) &
1255                    &                    +     e2t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) &
1256                    &                      ) / e2v(ji,jj)   +   e3v_0(ji,jj,jk)
1257                  END SELECT
1258               END DO
1259            END DO
1260         END DO
1261      ENDIF
1262         !                                             ! =====================
1263      IF( cp_cfg == "orca" .AND. jp_cfg == 05 ) THEN   ! ORCA R05 configuration
1264         !                                             ! =====================
1265         !
1266         ii0 = 563   ;   ii1 = 564        ! Gibraltar Strait (e2u was modified)
1267         ij0 = 327   ;   ij1 = 327
1268         DO jk = 1, jpkm1
1269            DO jj = mj0(ij0), mj1(ij1)
1270               DO ji = mi0(ii0), mi1(ii1)
1271                  SELECT CASE ( pout )
1272                  CASE( 'U' )
1273                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk)                                        &
1274                    &                    * (   e1t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) ) &
1275                    &                    +     e1t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) &
1276                    &                      ) / e1u(ji,jj)   +   e3u_0(ji,jj,jk)
1277                  CASE( 'F' )
1278                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk) * umask(ji,jj+1,jk)                    &
1279                    &                    * (   e1u(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3u_0(ji  ,jj,jk) ) &
1280                    &                    +     e1u(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3u_0(ji+1,jj,jk) ) &
1281                    &                      ) / e1f(ji,jj)   +   e3f_0(ji,jj,jk)
1282                  END SELECT
1283               END DO
1284            END DO
1285         END DO
1286         !
1287         ii0 = 627   ;   ii1 = 628        ! Bosphorus Strait (e2u was modified)
1288         ij0 = 343   ;   ij1 = 343
1289         DO jk = 1, jpkm1
1290            DO jj = mj0(ij0), mj1(ij1)
1291               DO ji = mi0(ii0), mi1(ii1)
1292                  SELECT CASE ( pout )
1293                  CASE( 'U' )
1294                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk)                                        & 
1295                    &                    * (   e1t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) ) &
1296                    &                    +     e1t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) &
1297                    &                      ) / e1u(ji,jj)   +   e3u_0(ji,jj,jk)
1298                  CASE( 'F' )
1299                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk) * umask(ji,jj+1,jk)                    & 
1300                    &                    * (   e1u(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3u_0(ji  ,jj,jk) ) &
1301                    &                    +     e1u(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3u_0(ji+1,jj,jk) ) &
1302                    &                      ) / e1f(ji,jj)   +   e3f_0(ji,jj,jk)
1303                  END SELECT
1304               END DO
1305            END DO
1306         END DO
1307         !
1308         ii0 =  93   ;   ii1 =  94        ! Sumba Strait (e2u was modified)
1309         ij0 = 232   ;   ij1 = 232
1310         DO jk = 1, jpkm1
1311            DO jj = mj0(ij0), mj1(ij1)
1312               DO ji = mi0(ii0), mi1(ii1)
1313                  SELECT CASE ( pout )
1314                  CASE( 'U' )
1315                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk)                                        &
1316                    &                    * (   e1t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) ) &
1317                    &                    +     e1t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) &
1318                    &                      ) / e1u(ji,jj)   +   e3u_0(ji,jj,jk)
1319                  CASE( 'F' )
1320                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk) * umask(ji,jj+1,jk)                    &
1321                    &                    * (   e1u(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3u_0(ji  ,jj,jk) ) &
1322                    &                    +     e1u(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3u_0(ji+1,jj,jk) ) &
1323                    &                      ) / e1f(ji,jj)   +   e3f_0(ji,jj,jk)
1324                  END SELECT
1325               END DO
1326            END DO
1327         END DO
1328         !
1329         ii0 = 103   ;   ii1 = 103        ! Ombai Strait (e2u was modified)
1330         ij0 = 232   ;   ij1 = 232
1331         DO jk = 1, jpkm1
1332            DO jj = mj0(ij0), mj1(ij1)
1333               DO ji = mi0(ii0), mi1(ii1)
1334                  SELECT CASE ( pout )
1335                  CASE( 'U' )
1336                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk)                                        &
1337                    &                    * (   e1t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) ) &
1338                    &                    +     e1t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) &
1339                    &                      ) / e1u(ji,jj)   +   e3u_0(ji,jj,jk)
1340                  CASE( 'F' )
1341                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk) * umask(ji,jj+1,jk)                    &
1342                    &                    * (   e1u(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3u_0(ji  ,jj,jk) ) &
1343                    &                    +     e1u(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3u_0(ji+1,jj,jk) ) &
1344                    &                      ) / e1f(ji,jj)   +   e3f_0(ji,jj,jk)
1345                  END SELECT
1346               END DO
1347            END DO
1348         END DO
1349         !
1350         ii0 =  15   ;   ii1 =  15        ! Palk Strait (e2u was modified)
1351         ij0 = 270   ;   ij1 = 270
1352         DO jk = 1, jpkm1
1353            DO jj = mj0(ij0), mj1(ij1)
1354               DO ji = mi0(ii0), mi1(ii1)
1355                  SELECT CASE ( pout )
1356                  CASE( 'U' )
1357                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk)                                        &
1358                    &                    * (   e1t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) ) &
1359                    &                    +     e1t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) &
1360                    &                      ) / e1u(ji,jj)   +   e3u_0(ji,jj,jk)
1361                  CASE( 'F' )
1362                     pe3_out(ji,jj,jk) = 0.5_wp * umask(ji,jj,jk) * umask(ji,jj+1,jk)                    &
1363                    &                    * (   e1u(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3u_0(ji  ,jj,jk) ) &
1364                    &                    +     e1u(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3u_0(ji+1,jj,jk) ) &
1365                    &                      ) / e1f(ji,jj)   +   e3f_0(ji,jj,jk)
1366                  END SELECT
1367               END DO
1368            END DO
1369         END DO
1370         !
1371         ii0 =  87   ;   ii1 =  87        ! Lombok Strait (e1v was modified)
1372         ij0 = 232   ;   ij1 = 233
1373         DO jk = 1, jpkm1
1374            DO jj = mj0(ij0), mj1(ij1)
1375               DO ji = mi0(ii0), mi1(ii1)
1376                  SELECT CASE ( pout )
1377                  CASE( 'V' )
1378                     pe3_out(ji,jj,jk) = 0.5_wp * vmask(ji,jj,jk)                                        &
1379                    &                    * (   e2t(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3t_0(ji,jj  ,jk) ) &
1380                    &                    +     e2t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) &
1381                    &                      ) / e2v(ji,jj)   +   e3v_0(ji,jj,jk)
1382                  END SELECT
1383               END DO
1384            END DO
1385         END DO
1386         !
1387         ii0 = 662   ;   ii1 = 662        ! Bab el Mandeb (e1v was modified)
1388         ij0 = 276   ;   ij1 = 276
1389         DO jk = 1, jpkm1
1390            DO jj = mj0(ij0), mj1(ij1)
1391               DO ji = mi0(ii0), mi1(ii1)
1392                  SELECT CASE ( pout )
1393                  CASE( 'V' )
1394                     pe3_out(ji,jj,jk) = 0.5_wp * vmask(ji,jj,jk)                                        &
1395                    &                    * (   e2t(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3t_0(ji,jj  ,jk) ) &
1396                    &                    +     e2t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) &
1397                    &                      ) / e2v(ji,jj)   +   e3v_0(ji,jj,jk)
1398                  END SELECT
1399               END DO
1400            END DO
1401         END DO
1402      ENDIF
1403   END SUBROUTINE dom_vvl_orca_fix
1404
1405   !!======================================================================
1406END MODULE domvvl
1407
1408
1409
Note: See TracBrowser for help on using the repository browser.