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

source: branches/UKMO/dev_r5518_Surge_Modelling/NEMOGCM/NEMO/OPA_SRC/DOM/domvvl.F90 @ 6510

Last change on this file since 6510 was 5954, checked in by rfurner, 9 years ago

added surge code from 2014_Surge_Modelling branch

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