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

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

source: NEMO/branches/2018/dev_r10164_HPC09_ESIWACE_PREP_MERGE/src/OCE/DOM/domvvl.F90 @ 10170

Last change on this file since 10170 was 10170, checked in by smasson, 5 years ago

dev_r10164_HPC09_ESIWACE_PREP_MERGE: action 2a: add report calls of lbc_lnk, see #2133

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