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_r7573_xios_write/NEMOGCM/NEMO/OPA_SRC/DOM – NEMO

source: branches/UKMO/dev_r7573_xios_write/NEMOGCM/NEMO/OPA_SRC/DOM/domvvl.F90 @ 8079

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

#1882 a first working version with XIOS writing restart file. Works with MO suite u-am389

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