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/trunk/src/OCE/DOM – NEMO

source: NEMO/trunk/src/OCE/DOM/domvvl.F90 @ 12108

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

trunk: merge dev_r10984_HPC-13 into the trunk

  • 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        )
[10425]81         CALL mpp_sum ( 'domvvl', dom_vvl_alloc )
82         IF( dom_vvl_alloc /= 0 )   CALL ctl_stop( 'STOP', 'dom_vvl_alloc: failed to allocate arrays' )
[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 )
[10425]88         CALL mpp_sum ( 'domvvl', dom_vvl_alloc )
89         IF( dom_vvl_alloc /= 0 )   CALL ctl_stop( 'STOP', 'dom_vvl_alloc: failed to allocate arrays' )
[4292]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
[10213]236            IF( cn_cfg == "orca" .OR. cn_cfg == "ORCA" ) THEN
237               IF( nn_cfg == 3 ) THEN   ! ORCA2: Suppress ztilde in the Foxe Basin for ORCA2
238                  ii0 = 103   ;   ii1 = 111       
239                  ij0 = 128   ;   ij1 = 135   ;   
240                  frq_rst_e3t( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) =  0.0_wp
241                  frq_rst_hdv( mi0(ii0):mi1(ii1) , mj0(ij0):mj1(ij1) ) =  1.e0_wp / rdt
242               ENDIF
[4292]243            ENDIF
244         ENDIF
245      ENDIF
[6140]246      !
[9367]247      IF(lwxios) THEN
248! define variables in restart file when writing with XIOS
249         CALL iom_set_rstw_var_active('e3t_b')
250         CALL iom_set_rstw_var_active('e3t_n')
251         !                                           ! ----------------------- !
252         IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN  ! z_tilde and layer cases !
253            !                                        ! ----------------------- !
254            CALL iom_set_rstw_var_active('tilde_e3t_b')
255            CALL iom_set_rstw_var_active('tilde_e3t_n')
256         END IF
257         !                                           ! -------------!   
258         IF( ln_vvl_ztilde ) THEN                    ! z_tilde case !
259            !                                        ! ------------ !
260            CALL iom_set_rstw_var_active('hdiv_lf')
261         ENDIF
262         !
263      ENDIF
264      !
[4292]265   END SUBROUTINE dom_vvl_init
266
267
[4338]268   SUBROUTINE dom_vvl_sf_nxt( kt, kcall ) 
[4292]269      !!----------------------------------------------------------------------
270      !!                ***  ROUTINE dom_vvl_sf_nxt  ***
271      !!                   
272      !! ** Purpose :  - compute the after scale factors used in tra_zdf, dynnxt,
273      !!                 tranxt and dynspg routines
274      !!
275      !! ** Method  :  - z_star case:  Repartition of ssh INCREMENT proportionnaly to the level thickness.
276      !!               - z_tilde_case: after scale factor increment =
277      !!                                    high frequency part of horizontal divergence
278      !!                                  + retsoring towards the background grid
279      !!                                  + thickness difusion
280      !!                               Then repartition of ssh INCREMENT proportionnaly
281      !!                               to the "baroclinic" level thickness.
282      !!
283      !! ** Action  :  - hdiv_lf    : restoring towards full baroclinic divergence in z_tilde case
284      !!               - tilde_e3t_a: after increment of vertical scale factor
285      !!                              in z_tilde case
[6140]286      !!               - e3(t/u/v)_a
[4292]287      !!
288      !! Reference  : Leclair, M., and Madec, G. 2011, Ocean Modelling.
289      !!----------------------------------------------------------------------
[6140]290      INTEGER, INTENT( in )           ::   kt      ! time step
291      INTEGER, INTENT( in ), OPTIONAL ::   kcall   ! optional argument indicating call sequence
292      !
293      INTEGER                ::   ji, jj, jk            ! dummy loop indices
294      INTEGER , DIMENSION(3) ::   ijk_max, ijk_min      ! temporary integers
295      REAL(wp)               ::   z2dt, z_tmin, z_tmax  ! local scalars
296      LOGICAL                ::   ll_do_bclinic         ! local logical
[9019]297      REAL(wp), DIMENSION(jpi,jpj)     ::   zht, z_scale, zwu, zwv, zhdiv
298      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   ze3t
[4292]299      !!----------------------------------------------------------------------
[6140]300      !
301      IF( ln_linssh )   RETURN      ! No calculation in linear free surface
302      !
[9019]303      IF( ln_timing )   CALL timing_start('dom_vvl_sf_nxt')
[6140]304      !
305      IF( kt == nit000 ) THEN
[4292]306         IF(lwp) WRITE(numout,*)
307         IF(lwp) WRITE(numout,*) 'dom_vvl_sf_nxt : compute after scale factors'
308         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~'
309      ENDIF
310
[4338]311      ll_do_bclinic = .TRUE.
312      IF( PRESENT(kcall) ) THEN
[6140]313         IF( kcall == 2 .AND. ln_vvl_ztilde )   ll_do_bclinic = .FALSE.
[4338]314      ENDIF
315
[4292]316      ! ******************************* !
317      ! After acale factors at t-points !
318      ! ******************************* !
[4338]319      !                                                ! --------------------------------------------- !
[6140]320      !                                                ! z_star coordinate and barotropic z-tilde part !
[4338]321      !                                                ! --------------------------------------------- !
[6140]322      !
[7753]323      z_scale(:,:) = ( ssha(:,:) - sshb(:,:) ) * ssmask(:,:) / ( ht_0(:,:) + sshn(:,:) + 1. - ssmask(:,:) )
[4338]324      DO jk = 1, jpkm1
[7753]325         ! formally this is the same as e3t_a = e3t_0*(1+ssha/ht_0)
326         e3t_a(:,:,jk) = e3t_b(:,:,jk) + e3t_n(:,:,jk) * z_scale(:,:) * tmask(:,:,jk)
[4338]327      END DO
[6140]328      !
[11415]329      IF( (ln_vvl_ztilde .OR. ln_vvl_layer) .AND. ll_do_bclinic ) THEN   ! z_tilde or layer coordinate !
330         !                                                               ! ------baroclinic part------ !
[4292]331         ! I - initialization
332         ! ==================
333
334         ! 1 - barotropic divergence
335         ! -------------------------
[7753]336         zhdiv(:,:) = 0._wp
337         zht(:,:)   = 0._wp
[4292]338         DO jk = 1, jpkm1
[7753]339            zhdiv(:,:) = zhdiv(:,:) + e3t_n(:,:,jk) * hdivn(:,:,jk)
340            zht  (:,:) = zht  (:,:) + e3t_n(:,:,jk) * tmask(:,:,jk)
[592]341         END DO
[7753]342         zhdiv(:,:) = zhdiv(:,:) / ( zht(:,:) + 1. - tmask_i(:,:) )
[2528]343
[4292]344         ! 2 - Low frequency baroclinic horizontal divergence  (z-tilde case only)
345         ! --------------------------------------------------
346         IF( ln_vvl_ztilde ) THEN
[6140]347            IF( kt > nit000 ) THEN
[4292]348               DO jk = 1, jpkm1
[7753]349                  hdiv_lf(:,:,jk) = hdiv_lf(:,:,jk) - rdt * frq_rst_hdv(:,:)   &
350                     &          * ( hdiv_lf(:,:,jk) - e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) )
[4292]351               END DO
352            ENDIF
[6140]353         ENDIF
[3294]354
[4292]355         ! II - after z_tilde increments of vertical scale factors
356         ! =======================================================
[7753]357         tilde_e3t_a(:,:,:) = 0._wp  ! tilde_e3t_a used to store tendency terms
[4292]358
359         ! 1 - High frequency divergence term
360         ! ----------------------------------
361         IF( ln_vvl_ztilde ) THEN     ! z_tilde case
362            DO jk = 1, jpkm1
[7753]363               tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - ( e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) - hdiv_lf(:,:,jk) )
[4292]364            END DO
365         ELSE                         ! layer case
366            DO jk = 1, jpkm1
[7753]367               tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) -   e3t_n(:,:,jk) * ( hdivn(:,:,jk) - zhdiv(:,:) ) * tmask(:,:,jk)
[4292]368            END DO
[6140]369         ENDIF
[4292]370
371         ! 2 - Restoring term (z-tilde case only)
372         ! ------------------
373         IF( ln_vvl_ztilde ) THEN
374            DO jk = 1, jpk
[7753]375               tilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - frq_rst_e3t(:,:) * tilde_e3t_b(:,:,jk)
[4292]376            END DO
[6140]377         ENDIF
[4292]378
379         ! 3 - Thickness diffusion term
380         ! ----------------------------
[7753]381         zwu(:,:) = 0._wp
382         zwv(:,:) = 0._wp
[6140]383         DO jk = 1, jpkm1        ! a - first derivative: diffusive fluxes
[4292]384            DO jj = 1, jpjm1
385               DO ji = 1, fs_jpim1   ! vector opt.
[5836]386                  un_td(ji,jj,jk) = rn_ahe3 * umask(ji,jj,jk) * e2_e1u(ji,jj)           &
387                     &            * ( tilde_e3t_b(ji,jj,jk) - tilde_e3t_b(ji+1,jj  ,jk) )
388                  vn_td(ji,jj,jk) = rn_ahe3 * vmask(ji,jj,jk) * e1_e2v(ji,jj)           & 
389                     &            * ( tilde_e3t_b(ji,jj,jk) - tilde_e3t_b(ji  ,jj+1,jk) )
[4292]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
[6140]395         DO jj = 1, jpj          ! b - correction for last oceanic u-v points
[4292]396            DO ji = 1, jpi
397               un_td(ji,jj,mbku(ji,jj)) = un_td(ji,jj,mbku(ji,jj)) - zwu(ji,jj)
398               vn_td(ji,jj,mbkv(ji,jj)) = vn_td(ji,jj,mbkv(ji,jj)) - zwv(ji,jj)
399            END DO
400         END DO
[6140]401         DO jk = 1, jpkm1        ! c - second derivative: divergence of diffusive fluxes
[4292]402            DO jj = 2, jpjm1
403               DO ji = fs_2, fs_jpim1   ! vector opt.
404                  tilde_e3t_a(ji,jj,jk) = tilde_e3t_a(ji,jj,jk) + (   un_td(ji-1,jj  ,jk) - un_td(ji,jj,jk)    &
405                     &                                          +     vn_td(ji  ,jj-1,jk) - vn_td(ji,jj,jk)    &
[5836]406                     &                                            ) * r1_e1e2t(ji,jj)
[4292]407               END DO
408            END DO
409         END DO
[6140]410         !                       ! d - thickness diffusion transport: boundary conditions
411         !                             (stored for tracer advction and continuity equation)
[10425]412         CALL lbc_lnk_multi( 'domvvl', un_td , 'U' , -1._wp, vn_td , 'V' , -1._wp)
[4292]413
414         ! 4 - Time stepping of baroclinic scale factors
415         ! ---------------------------------------------
416         ! Leapfrog time stepping
417         ! ~~~~~~~~~~~~~~~~~~~~~~
418         IF( neuler == 0 .AND. kt == nit000 ) THEN
419            z2dt =  rdt
420         ELSE
421            z2dt = 2.0_wp * rdt
422         ENDIF
[10425]423         CALL lbc_lnk( 'domvvl', tilde_e3t_a(:,:,:), 'T', 1._wp )
[7753]424         tilde_e3t_a(:,:,:) = tilde_e3t_b(:,:,:) + z2dt * tmask(:,:,:) * tilde_e3t_a(:,:,:)
[4292]425
426         ! Maximum deformation control
427         ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~
[7753]428         ze3t(:,:,jpk) = 0._wp
[4292]429         DO jk = 1, jpkm1
[7753]430            ze3t(:,:,jk) = tilde_e3t_a(:,:,jk) / e3t_0(:,:,jk) * tmask(:,:,jk) * tmask_i(:,:)
[4292]431         END DO
432         z_tmax = MAXVAL( ze3t(:,:,:) )
[10425]433         CALL mpp_max( 'domvvl', z_tmax )                 ! max over the global domain
[4292]434         z_tmin = MINVAL( ze3t(:,:,:) )
[10425]435         CALL mpp_min( 'domvvl', z_tmin )                 ! min over the global domain
[4292]436         ! - ML - test: for the moment, stop simulation for too large e3_t variations
[6140]437         IF( ( z_tmax >  rn_zdef_max ) .OR. ( z_tmin < - rn_zdef_max ) ) THEN
[4292]438            IF( lk_mpp ) THEN
[10425]439               CALL mpp_maxloc( 'domvvl', ze3t, tmask, z_tmax, ijk_max )
440               CALL mpp_minloc( 'domvvl', ze3t, tmask, z_tmin, ijk_min )
[4292]441            ELSE
442               ijk_max = MAXLOC( ze3t(:,:,:) )
443               ijk_max(1) = ijk_max(1) + nimpp - 1
444               ijk_max(2) = ijk_max(2) + njmpp - 1
445               ijk_min = MINLOC( ze3t(:,:,:) )
446               ijk_min(1) = ijk_min(1) + nimpp - 1
447               ijk_min(2) = ijk_min(2) + njmpp - 1
448            ENDIF
449            IF (lwp) THEN
450               WRITE(numout, *) 'MAX( tilde_e3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmax
451               WRITE(numout, *) 'at i, j, k=', ijk_max
452               WRITE(numout, *) 'MIN( tilde_e3t_a(:,:,:) / e3t_0(:,:,:) ) =', z_tmin
453               WRITE(numout, *) 'at i, j, k=', ijk_min           
[10425]454               CALL ctl_stop( 'STOP', 'MAX( ABS( tilde_e3t_a(:,:,: ) ) / e3t_0(:,:,:) ) too high')
[4292]455            ENDIF
456         ENDIF
457         ! - ML - end test
458         ! - ML - Imposing these limits will cause a baroclinicity error which is corrected for below
[7753]459         tilde_e3t_a(:,:,:) = MIN( tilde_e3t_a(:,:,:),   rn_zdef_max * e3t_0(:,:,:) )
460         tilde_e3t_a(:,:,:) = MAX( tilde_e3t_a(:,:,:), - rn_zdef_max * e3t_0(:,:,:) )
[4292]461
[4338]462         !
463         ! "tilda" change in the after scale factor
[4292]464         ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
[4338]465         DO jk = 1, jpkm1
[7753]466            dtilde_e3t_a(:,:,jk) = tilde_e3t_a(:,:,jk) - tilde_e3t_b(:,:,jk)
[4338]467         END DO
[4292]468         ! III - Barotropic repartition of the sea surface height over the baroclinic profile
469         ! ==================================================================================
[4338]470         ! add ( ssh increment + "baroclinicity error" ) proportionly to e3t(n)
[4292]471         ! - ML - baroclinicity error should be better treated in the future
472         !        i.e. locally and not spread over the water column.
473         !        (keep in mind that the idea is to reduce Eulerian velocity as much as possible)
[7753]474         zht(:,:) = 0.
[4292]475         DO jk = 1, jpkm1
[7753]476            zht(:,:)  = zht(:,:) + tilde_e3t_a(:,:,jk) * tmask(:,:,jk)
[4292]477         END DO
[7753]478         z_scale(:,:) =  - zht(:,:) / ( ht_0(:,:) + sshn(:,:) + 1. - ssmask(:,:) )
[4292]479         DO jk = 1, jpkm1
[7753]480            dtilde_e3t_a(:,:,jk) = dtilde_e3t_a(:,:,jk) + e3t_n(:,:,jk) * z_scale(:,:) * tmask(:,:,jk)
[4292]481         END DO
[7753]482
[4292]483      ENDIF
484
[4338]485      IF( ln_vvl_ztilde .OR. ln_vvl_layer )  THEN   ! z_tilde or layer coordinate !
486      !                                           ! ---baroclinic part--------- !
487         DO jk = 1, jpkm1
[7753]488            e3t_a(:,:,jk) = e3t_a(:,:,jk) + dtilde_e3t_a(:,:,jk) * tmask(:,:,jk)
[4338]489         END DO
490      ENDIF
491
492      IF( ln_vvl_dbg .AND. .NOT. ll_do_bclinic ) THEN   ! - ML - test: control prints for debuging
[4292]493         !
494         IF( lwp ) WRITE(numout, *) 'kt =', kt
495         IF ( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
496            z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( zht(:,:) ) )
[10425]497            CALL mpp_max( 'domvvl', z_tmax )                             ! max over the global domain
[4292]498            IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(SUM(tilde_e3t_a))) =', z_tmax
499         END IF
500         !
[7753]501         zht(:,:) = 0.0_wp
[4292]502         DO jk = 1, jpkm1
[7753]503            zht(:,:) = zht(:,:) + e3t_n(:,:,jk) * tmask(:,:,jk)
[4292]504         END DO
505         z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + sshn(:,:) - zht(:,:) ) )
[10425]506         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
[6140]507         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+sshn-SUM(e3t_n))) =', z_tmax
[4292]508         !
[7753]509         zht(:,:) = 0.0_wp
[4292]510         DO jk = 1, jpkm1
[7753]511            zht(:,:) = zht(:,:) + e3t_a(:,:,jk) * tmask(:,:,jk)
[4292]512         END DO
513         z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + ssha(:,:) - zht(:,:) ) )
[10425]514         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
[6140]515         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+ssha-SUM(e3t_a))) =', z_tmax
[4292]516         !
[7753]517         zht(:,:) = 0.0_wp
[4292]518         DO jk = 1, jpkm1
[7753]519            zht(:,:) = zht(:,:) + e3t_b(:,:,jk) * tmask(:,:,jk)
[4292]520         END DO
521         z_tmax = MAXVAL( tmask(:,:,1) * tmask_i(:,:) * ABS( ht_0(:,:) + sshb(:,:) - zht(:,:) ) )
[10425]522         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
[6140]523         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ht_0+sshb-SUM(e3t_b))) =', z_tmax
[4292]524         !
525         z_tmax = MAXVAL( tmask(:,:,1) *  ABS( sshb(:,:) ) )
[10425]526         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
[4292]527         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(sshb))) =', z_tmax
528         !
529         z_tmax = MAXVAL( tmask(:,:,1) *  ABS( sshn(:,:) ) )
[10425]530         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
[4292]531         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(sshn))) =', z_tmax
532         !
533         z_tmax = MAXVAL( tmask(:,:,1) *  ABS( ssha(:,:) ) )
[10425]534         CALL mpp_max( 'domvvl', z_tmax )                                ! max over the global domain
[4292]535         IF( lwp    ) WRITE(numout, *) kt,' MAXVAL(abs(ssha))) =', z_tmax
536      END IF
537
538      ! *********************************** !
539      ! After scale factors at u- v- points !
540      ! *********************************** !
541
[6140]542      CALL dom_vvl_interpol( e3t_a(:,:,:), e3u_a(:,:,:), 'U' )
543      CALL dom_vvl_interpol( e3t_a(:,:,:), e3v_a(:,:,:), 'V' )
[4292]544
[4370]545      ! *********************************** !
546      ! After depths at u- v points         !
547      ! *********************************** !
548
[7753]549      hu_a(:,:) = e3u_a(:,:,1) * umask(:,:,1)
550      hv_a(:,:) = e3v_a(:,:,1) * vmask(:,:,1)
[6140]551      DO jk = 2, jpkm1
[7753]552         hu_a(:,:) = hu_a(:,:) + e3u_a(:,:,jk) * umask(:,:,jk)
553         hv_a(:,:) = hv_a(:,:) + e3v_a(:,:,jk) * vmask(:,:,jk)
[4370]554      END DO
555      !                                        ! Inverse of the local depth
[6140]556!!gm BUG ?  don't understand the use of umask_i here .....
[7753]557      r1_hu_a(:,:) = ssumask(:,:) / ( hu_a(:,:) + 1._wp - ssumask(:,:) )
558      r1_hv_a(:,:) = ssvmask(:,:) / ( hv_a(:,:) + 1._wp - ssvmask(:,:) )
[6140]559      !
[9019]560      IF( ln_timing )   CALL timing_stop('dom_vvl_sf_nxt')
[6140]561      !
[4292]562   END SUBROUTINE dom_vvl_sf_nxt
563
564
565   SUBROUTINE dom_vvl_sf_swp( kt )
[3294]566      !!----------------------------------------------------------------------
[4292]567      !!                ***  ROUTINE dom_vvl_sf_swp  ***
[3294]568      !!                   
[4292]569      !! ** Purpose :  compute time filter and swap of scale factors
570      !!               compute all depths and related variables for next time step
571      !!               write outputs and restart file
[3294]572      !!
[4292]573      !! ** Method  :  - swap of e3t with trick for volume/tracer conservation
574      !!               - reconstruct scale factor at other grid points (interpolate)
575      !!               - recompute depths and water height fields
576      !!
[6140]577      !! ** Action  :  - e3t_(b/n), tilde_e3t_(b/n) and e3(u/v)_n ready for next time step
[4292]578      !!               - Recompute:
[6140]579      !!                    e3(u/v)_b       
580      !!                    e3w_n           
581      !!                    e3(u/v)w_b     
582      !!                    e3(u/v)w_n     
583      !!                    gdept_n, gdepw_n  and gde3w_n
[4292]584      !!                    h(u/v) and h(u/v)r
585      !!
586      !! Reference  : Leclair, M., and G. Madec, 2009, Ocean Modelling.
587      !!              Leclair, M., and G. Madec, 2011, Ocean Modelling.
[3294]588      !!----------------------------------------------------------------------
[6140]589      INTEGER, INTENT( in ) ::   kt   ! time step
590      !
591      INTEGER  ::   ji, jj, jk   ! dummy loop indices
592      REAL(wp) ::   zcoef        ! local scalar
[3294]593      !!----------------------------------------------------------------------
[6140]594      !
595      IF( ln_linssh )   RETURN      ! No calculation in linear free surface
596      !
[9019]597      IF( ln_timing )   CALL timing_start('dom_vvl_sf_swp')
[3294]598      !
[4292]599      IF( kt == nit000 )   THEN
600         IF(lwp) WRITE(numout,*)
601         IF(lwp) WRITE(numout,*) 'dom_vvl_sf_swp : - time filter and swap of scale factors'
602         IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~~   - interpolate scale factors and compute depths for next time step'
[3294]603      ENDIF
[4292]604      !
605      ! Time filter and swap of scale factors
606      ! =====================================
[6140]607      ! - ML - e3(t/u/v)_b are allready computed in dynnxt.
[4292]608      IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN
609         IF( neuler == 0 .AND. kt == nit000 ) THEN
[7753]610            tilde_e3t_b(:,:,:) = tilde_e3t_n(:,:,:)
[4292]611         ELSE
[7753]612            tilde_e3t_b(:,:,:) = tilde_e3t_n(:,:,:) & 
613            &         + atfp * ( tilde_e3t_b(:,:,:) - 2.0_wp * tilde_e3t_n(:,:,:) + tilde_e3t_a(:,:,:) )
[4292]614         ENDIF
[7753]615         tilde_e3t_n(:,:,:) = tilde_e3t_a(:,:,:)
[4292]616      ENDIF
[7753]617      gdept_b(:,:,:) = gdept_n(:,:,:)
618      gdepw_b(:,:,:) = gdepw_n(:,:,:)
[4488]619
[7753]620      e3t_n(:,:,:) = e3t_a(:,:,:)
621      e3u_n(:,:,:) = e3u_a(:,:,:)
622      e3v_n(:,:,:) = e3v_a(:,:,:)
623
[4292]624      ! Compute all missing vertical scale factor and depths
625      ! ====================================================
626      ! Horizontal scale factor interpolations
627      ! --------------------------------------
[6140]628      ! - ML - e3u_b and e3v_b are allready computed in dynnxt
[4370]629      ! - JC - hu_b, hv_b, hur_b, hvr_b also
[6140]630     
631      CALL dom_vvl_interpol( e3u_n(:,:,:), e3f_n(:,:,:), 'F'  )
632     
[4292]633      ! Vertical scale factor interpolations
[6140]634      CALL dom_vvl_interpol( e3t_n(:,:,:),  e3w_n(:,:,:), 'W'  )
635      CALL dom_vvl_interpol( e3u_n(:,:,:), e3uw_n(:,:,:), 'UW' )
636      CALL dom_vvl_interpol( e3v_n(:,:,:), e3vw_n(:,:,:), 'VW' )
637      CALL dom_vvl_interpol( e3t_b(:,:,:),  e3w_b(:,:,:), 'W'  )
638      CALL dom_vvl_interpol( e3u_b(:,:,:), e3uw_b(:,:,:), 'UW' )
639      CALL dom_vvl_interpol( e3v_b(:,:,:), e3vw_b(:,:,:), 'VW' )
[5120]640
[6140]641      ! t- and w- points depth (set the isf depth as it is in the initial step)
[7753]642      gdept_n(:,:,1) = 0.5_wp * e3w_n(:,:,1)
643      gdepw_n(:,:,1) = 0.0_wp
644      gde3w_n(:,:,1) = gdept_n(:,:,1) - sshn(:,:)
[5120]645      DO jk = 2, jpk
646         DO jj = 1,jpj
647            DO ji = 1,jpi
648              !    zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))   ! 0 everywhere tmask = wmask, ie everywhere expect at jk = mikt
649                                                                 ! 1 for jk = mikt
650               zcoef = (tmask(ji,jj,jk) - wmask(ji,jj,jk))
[6140]651               gdepw_n(ji,jj,jk) = gdepw_n(ji,jj,jk-1) + e3t_n(ji,jj,jk-1)
652               gdept_n(ji,jj,jk) =    zcoef  * ( gdepw_n(ji,jj,jk  ) + 0.5 * e3w_n(ji,jj,jk) )  &
653                   &             + (1-zcoef) * ( gdept_n(ji,jj,jk-1) +       e3w_n(ji,jj,jk) ) 
654               gde3w_n(ji,jj,jk) = gdept_n(ji,jj,jk) - sshn(ji,jj)
[4990]655            END DO
656         END DO
[4292]657      END DO
[5120]658
[6140]659      ! Local depth and Inverse of the local depth of the water
660      ! -------------------------------------------------------
[7753]661      hu_n(:,:) = hu_a(:,:)   ;   r1_hu_n(:,:) = r1_hu_a(:,:)
662      hv_n(:,:) = hv_a(:,:)   ;   r1_hv_n(:,:) = r1_hv_a(:,:)
663      !
664      ht_n(:,:) = e3t_n(:,:,1) * tmask(:,:,1)
[6140]665      DO jk = 2, jpkm1
[7753]666         ht_n(:,:) = ht_n(:,:) + e3t_n(:,:,jk) * tmask(:,:,jk)
[4370]667      END DO
[7753]668
[4292]669      ! write restart file
670      ! ==================
[9019]671      IF( lrst_oce  )   CALL dom_vvl_rst( kt, 'WRITE' )
[4292]672      !
[9019]673      IF( ln_timing )   CALL timing_stop('dom_vvl_sf_swp')
[6140]674      !
[4292]675   END SUBROUTINE dom_vvl_sf_swp
676
677
678   SUBROUTINE dom_vvl_interpol( pe3_in, pe3_out, pout )
679      !!---------------------------------------------------------------------
680      !!                  ***  ROUTINE dom_vvl__interpol  ***
681      !!
682      !! ** Purpose :   interpolate scale factors from one grid point to another
683      !!
684      !! ** Method  :   e3_out = e3_0 + interpolation(e3_in - e3_0)
685      !!                - horizontal interpolation: grid cell surface averaging
686      !!                - vertical interpolation: simple averaging
687      !!----------------------------------------------------------------------
[5836]688      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) ::  pe3_in    ! input e3 to be interpolated
689      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) ::  pe3_out   ! output interpolated e3
690      CHARACTER(LEN=*)                , INTENT(in   ) ::  pout      ! grid point of out scale factors
691      !                                                             !   =  'U', 'V', 'W, 'F', 'UW' or 'VW'
692      !
[6152]693      INTEGER ::   ji, jj, jk                                       ! dummy loop indices
[9023]694      REAL(wp) ::  zlnwd                                            ! =1./0. when ln_wd_il = T/F
[4292]695      !!----------------------------------------------------------------------
[5836]696      !
[9023]697      IF(ln_wd_il) THEN
[6152]698        zlnwd = 1.0_wp
699      ELSE
700        zlnwd = 0.0_wp
701      END IF
702      !
[5836]703      SELECT CASE ( pout )    !==  type of interpolation  ==!
[4292]704         !
[5836]705      CASE( 'U' )                   !* from T- to U-point : hor. surface weighted mean
[4292]706         DO jk = 1, jpk
707            DO jj = 1, jpjm1
708               DO ji = 1, fs_jpim1   ! vector opt.
[6152]709                  pe3_out(ji,jj,jk) = 0.5_wp * (  umask(ji,jj,jk) * (1.0_wp - zlnwd) + zlnwd ) * r1_e1e2u(ji,jj)   &
[5836]710                     &                       * (   e1e2t(ji  ,jj) * ( pe3_in(ji  ,jj,jk) - e3t_0(ji  ,jj,jk) )     &
711                     &                           + e1e2t(ji+1,jj) * ( pe3_in(ji+1,jj,jk) - e3t_0(ji+1,jj,jk) ) )
[4292]712               END DO
[2528]713            END DO
714         END DO
[10425]715         CALL lbc_lnk( 'domvvl', pe3_out(:,:,:), 'U', 1._wp )
[7753]716         pe3_out(:,:,:) = pe3_out(:,:,:) + e3u_0(:,:,:)
[5836]717         !
718      CASE( 'V' )                   !* from T- to V-point : hor. surface weighted mean
[4292]719         DO jk = 1, jpk
720            DO jj = 1, jpjm1
721               DO ji = 1, fs_jpim1   ! vector opt.
[6152]722                  pe3_out(ji,jj,jk) = 0.5_wp * ( vmask(ji,jj,jk)  * (1.0_wp - zlnwd) + zlnwd ) * r1_e1e2v(ji,jj)   &
[5836]723                     &                       * (   e1e2t(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3t_0(ji,jj  ,jk) )     &
724                     &                           + e1e2t(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3t_0(ji,jj+1,jk) ) )
[4292]725               END DO
726            END DO
727         END DO
[10425]728         CALL lbc_lnk( 'domvvl', pe3_out(:,:,:), 'V', 1._wp )
[7753]729         pe3_out(:,:,:) = pe3_out(:,:,:) + e3v_0(:,:,:)
[5836]730         !
731      CASE( 'F' )                   !* from U-point to F-point : hor. surface weighted mean
[4292]732         DO jk = 1, jpk
733            DO jj = 1, jpjm1
734               DO ji = 1, fs_jpim1   ! vector opt.
[6152]735                  pe3_out(ji,jj,jk) = 0.5_wp * (  umask(ji,jj,jk) * umask(ji,jj+1,jk) * (1.0_wp - zlnwd) + zlnwd ) &
736                     &                       *    r1_e1e2f(ji,jj)                                                  &
[5836]737                     &                       * (   e1e2u(ji,jj  ) * ( pe3_in(ji,jj  ,jk) - e3u_0(ji,jj  ,jk) )     &
738                     &                           + e1e2u(ji,jj+1) * ( pe3_in(ji,jj+1,jk) - e3u_0(ji,jj+1,jk) ) )
[4292]739               END DO
740            END DO
741         END DO
[10425]742         CALL lbc_lnk( 'domvvl', pe3_out(:,:,:), 'F', 1._wp )
[7753]743         pe3_out(:,:,:) = pe3_out(:,:,:) + e3f_0(:,:,:)
[5836]744         !
745      CASE( 'W' )                   !* from T- to W-point : vertical simple mean
746         !
[7753]747         pe3_out(:,:,1) = e3w_0(:,:,1) + pe3_in(:,:,1) - e3t_0(:,:,1)
[5836]748         ! - ML - The use of mask in this formulea enables the special treatment of the last w-point without indirect adressing
749!!gm BUG? use here wmask in case of ISF ?  to be checked
[4292]750         DO jk = 2, jpk
[7753]751            pe3_out(:,:,jk) = e3w_0(:,:,jk) + ( 1.0_wp - 0.5_wp * ( tmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) )   &
752               &                            * ( pe3_in(:,:,jk-1) - e3t_0(:,:,jk-1) )                               &
753               &                            +            0.5_wp * ( tmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd )     &
754               &                            * ( pe3_in(:,:,jk  ) - e3t_0(:,:,jk  ) )
[4292]755         END DO
[5836]756         !
757      CASE( 'UW' )                  !* from U- to UW-point : vertical simple mean
758         !
[7753]759         pe3_out(:,:,1) = e3uw_0(:,:,1) + pe3_in(:,:,1) - e3u_0(:,:,1)
[4292]760         ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without indirect adressing
[5836]761!!gm BUG? use here wumask in case of ISF ?  to be checked
[4292]762         DO jk = 2, jpk
[7753]763            pe3_out(:,:,jk) = e3uw_0(:,:,jk) + ( 1.0_wp - 0.5_wp * ( umask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) )  &
764               &                             * ( pe3_in(:,:,jk-1) - e3u_0(:,:,jk-1) )                              &
765               &                             +            0.5_wp * ( umask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd )    &
766               &                             * ( pe3_in(:,:,jk  ) - e3u_0(:,:,jk  ) )
[4292]767         END DO
[5836]768         !
769      CASE( 'VW' )                  !* from V- to VW-point : vertical simple mean
770         !
[7753]771         pe3_out(:,:,1) = e3vw_0(:,:,1) + pe3_in(:,:,1) - e3v_0(:,:,1)
[4292]772         ! - ML - The use of mask in this formaula enables the special treatment of the last w- point without indirect adressing
[5836]773!!gm BUG? use here wvmask in case of ISF ?  to be checked
[4292]774         DO jk = 2, jpk
[7753]775            pe3_out(:,:,jk) = e3vw_0(:,:,jk) + ( 1.0_wp - 0.5_wp * ( vmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd ) )  &
776               &                             * ( pe3_in(:,:,jk-1) - e3v_0(:,:,jk-1) )                              &
777               &                             +            0.5_wp * ( vmask(:,:,jk) * (1.0_wp - zlnwd) + zlnwd )    &
778               &                             * ( pe3_in(:,:,jk  ) - e3v_0(:,:,jk  ) )
[4292]779         END DO
780      END SELECT
781      !
782   END SUBROUTINE dom_vvl_interpol
783
[5836]784
[4292]785   SUBROUTINE dom_vvl_rst( kt, cdrw )
786      !!---------------------------------------------------------------------
787      !!                   ***  ROUTINE dom_vvl_rst  ***
788      !!                     
789      !! ** Purpose :   Read or write VVL file in restart file
790      !!
791      !! ** Method  :   use of IOM library
792      !!                if the restart does not contain vertical scale factors,
793      !!                they are set to the _0 values
794      !!                if the restart does not contain vertical scale factors increments (z_tilde),
795      !!                they are set to 0.
796      !!----------------------------------------------------------------------
797      INTEGER         , INTENT(in) ::   kt     ! ocean time-step
798      CHARACTER(len=*), INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
[5836]799      !
[6152]800      INTEGER ::   ji, jj, jk
[4292]801      INTEGER ::   id1, id2, id3, id4, id5     ! local integers
802      !!----------------------------------------------------------------------
803      !
804      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialise
805         !                                   ! ===============
806         IF( ln_rstart ) THEN                   !* Read the restart file
807            CALL rst_read_open                  !  open the restart file if necessary
[9367]808            CALL iom_get( numror, jpdom_autoglo, 'sshn'   , sshn, ldxios = lrxios    )
[4366]809            !
[6140]810            id1 = iom_varid( numror, 'e3t_b', ldstop = .FALSE. )
811            id2 = iom_varid( numror, 'e3t_n', ldstop = .FALSE. )
[4292]812            id3 = iom_varid( numror, 'tilde_e3t_b', ldstop = .FALSE. )
813            id4 = iom_varid( numror, 'tilde_e3t_n', ldstop = .FALSE. )
[4795]814            id5 = iom_varid( numror, 'hdiv_lf', ldstop = .FALSE. )
[4292]815            !                             ! --------- !
816            !                             ! all cases !
817            !                             ! --------- !
818            IF( MIN( id1, id2 ) > 0 ) THEN       ! all required arrays exist
[9367]819               CALL iom_get( numror, jpdom_autoglo, 'e3t_b', e3t_b(:,:,:), ldxios = lrxios )
820               CALL iom_get( numror, jpdom_autoglo, 'e3t_n', e3t_n(:,:,:), ldxios = lrxios )
[4990]821               ! needed to restart if land processor not computed
[6140]822               IF(lwp) write(numout,*) 'dom_vvl_rst : e3t_b and e3t_n found in restart files'
[4990]823               WHERE ( tmask(:,:,:) == 0.0_wp ) 
[6140]824                  e3t_n(:,:,:) = e3t_0(:,:,:)
825                  e3t_b(:,:,:) = e3t_0(:,:,:)
[4990]826               END WHERE
[4292]827               IF( neuler == 0 ) THEN
[6140]828                  e3t_b(:,:,:) = e3t_n(:,:,:)
[4292]829               ENDIF
830            ELSE IF( id1 > 0 ) THEN
[6140]831               IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : e3t_n not found in restart files'
832               IF(lwp) write(numout,*) 'e3t_n set equal to e3t_b.'
[4990]833               IF(lwp) write(numout,*) 'neuler is forced to 0'
[9367]834               CALL iom_get( numror, jpdom_autoglo, 'e3t_b', e3t_b(:,:,:), ldxios = lrxios )
[6140]835               e3t_n(:,:,:) = e3t_b(:,:,:)
[4990]836               neuler = 0
837            ELSE IF( id2 > 0 ) THEN
[6140]838               IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : e3t_b not found in restart files'
839               IF(lwp) write(numout,*) 'e3t_b set equal to e3t_n.'
[4490]840               IF(lwp) write(numout,*) 'neuler is forced to 0'
[9367]841               CALL iom_get( numror, jpdom_autoglo, 'e3t_n', e3t_n(:,:,:), ldxios = lrxios )
[6140]842               e3t_b(:,:,:) = e3t_n(:,:,:)
[4490]843               neuler = 0
844            ELSE
[6140]845               IF(lwp) write(numout,*) 'dom_vvl_rst WARNING : e3t_n not found in restart file'
[4490]846               IF(lwp) write(numout,*) 'Compute scale factor from sshn'
847               IF(lwp) write(numout,*) 'neuler is forced to 0'
[6140]848               DO jk = 1, jpk
849                  e3t_n(:,:,jk) =  e3t_0(:,:,jk) * ( ht_0(:,:) + sshn(:,:) ) &
850                      &                          / ( ht_0(:,:) + 1._wp - ssmask(:,:) ) * tmask(:,:,jk)   &
851                      &          + e3t_0(:,:,jk)                               * (1._wp -tmask(:,:,jk))
[4490]852               END DO
[6140]853               e3t_b(:,:,:) = e3t_n(:,:,:)
[4490]854               neuler = 0
[4292]855            ENDIF
856            !                             ! ----------- !
857            IF( ln_vvl_zstar ) THEN       ! z_star case !
858               !                          ! ----------- !
859               IF( MIN( id3, id4 ) > 0 ) THEN
860                  CALL ctl_stop( 'dom_vvl_rst: z_star cannot restart from a z_tilde or layer run' )
861               ENDIF
862               !                          ! ----------------------- !
863            ELSE                          ! z_tilde and layer cases !
864               !                          ! ----------------------- !
865               IF( MIN( id3, id4 ) > 0 ) THEN  ! all required arrays exist
[9367]866                  CALL iom_get( numror, jpdom_autoglo, 'tilde_e3t_b', tilde_e3t_b(:,:,:), ldxios = lrxios )
867                  CALL iom_get( numror, jpdom_autoglo, 'tilde_e3t_n', tilde_e3t_n(:,:,:), ldxios = lrxios )
[4292]868               ELSE                            ! one at least array is missing
869                  tilde_e3t_b(:,:,:) = 0.0_wp
870                  tilde_e3t_n(:,:,:) = 0.0_wp
871               ENDIF
872               !                          ! ------------ !
873               IF( ln_vvl_ztilde ) THEN   ! z_tilde case !
874                  !                       ! ------------ !
875                  IF( id5 > 0 ) THEN  ! required array exists
[9367]876                     CALL iom_get( numror, jpdom_autoglo, 'hdiv_lf', hdiv_lf(:,:,:), ldxios = lrxios )
[4292]877                  ELSE                ! array is missing
878                     hdiv_lf(:,:,:) = 0.0_wp
879                  ENDIF
880               ENDIF
881            ENDIF
882            !
883         ELSE                                   !* Initialize at "rest"
[7646]884            !
[9023]885
886            IF( ll_wd ) THEN   ! MJB ll_wd edits start here - these are essential
887               !
888               IF( cn_cfg == 'wad' ) THEN
889                  ! Wetting and drying test case
890                  CALL usr_def_istate( gdept_b, tmask, tsb, ub, vb, sshb  )
891                  tsn  (:,:,:,:) = tsb (:,:,:,:)       ! set now values from to before ones
892                  sshn (:,:)     = sshb(:,:)
893                  un   (:,:,:)   = ub  (:,:,:)
894                  vn   (:,:,:)   = vb  (:,:,:)
895               ELSE
896                  ! if not test case
897                  sshn(:,:) = -ssh_ref
898                  sshb(:,:) = -ssh_ref
899
900                  DO jj = 1, jpj
901                     DO ji = 1, jpi
902                        IF( ht_0(ji,jj)-ssh_ref <  rn_wdmin1 ) THEN ! if total depth is less than min depth
903
904                           sshb(ji,jj) = rn_wdmin1 - (ht_0(ji,jj) )
905                           sshn(ji,jj) = rn_wdmin1 - (ht_0(ji,jj) )
906                           ssha(ji,jj) = rn_wdmin1 - (ht_0(ji,jj) )
907                        ENDIF
908                     ENDDO
909                  ENDDO
910               ENDIF !If test case else
911
912               ! Adjust vertical metrics for all wad
[7646]913               DO jk = 1, jpk
[9023]914                  e3t_n(:,:,jk) =  e3t_0(:,:,jk) * ( ht_0(:,:) + sshn(:,:)  ) &
[7646]915                    &                            / ( ht_0(:,:) + 1._wp - ssmask(:,:) ) * tmask(:,:,jk)   &
[9023]916                    &            + e3t_0(:,:,jk) * ( 1._wp - tmask(:,:,jk) )
[7646]917               END DO
918               e3t_b(:,:,:) = e3t_n(:,:,:)
[9023]919
920               DO ji = 1, jpi
921                  DO jj = 1, jpj
922                     IF ( ht_0(ji,jj) .LE. 0.0 .AND. NINT( ssmask(ji,jj) ) .EQ. 1) THEN
923                       CALL ctl_stop( 'dom_vvl_rst: ht_0 must be positive at potentially wet points' )
924                     ENDIF
925                  END DO
926               END DO 
[7646]927               !
928            ELSE
929               !
[9059]930               ! Just to read set ssh in fact, called latter once vertical grid
931               ! is set up:
[9065]932!               CALL usr_def_istate( gdept_0, tmask, tsb, ub, vb, sshb  )
933!               !
934!               DO jk=1,jpk
935!                  e3t_b(:,:,jk) =  e3t_0(:,:,jk) * ( ht_0(:,:) + sshb(:,:) ) &
936!                     &            / ( ht_0(:,:) + 1._wp -ssmask(:,:) ) * tmask(:,:,jk)
937!               END DO
938!               e3t_n(:,:,:) = e3t_b(:,:,:)
939                sshn(:,:)=0._wp
940                e3t_n(:,:,:)=e3t_0(:,:,:)
941                e3t_b(:,:,:)=e3t_0(:,:,:)
[7646]942               !
[9023]943            END IF           ! end of ll_wd edits
[6152]944
[4292]945            IF( ln_vvl_ztilde .OR. ln_vvl_layer) THEN
[7646]946               tilde_e3t_b(:,:,:) = 0._wp
947               tilde_e3t_n(:,:,:) = 0._wp
948               IF( ln_vvl_ztilde ) hdiv_lf(:,:,:) = 0._wp
[4292]949            END IF
950         ENDIF
[5836]951         !
[4292]952      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
953         !                                   ! ===================
954         IF(lwp) WRITE(numout,*) '---- dom_vvl_rst ----'
[9367]955         IF( lwxios ) CALL iom_swap(      cwxios_context          )
[4292]956         !                                           ! --------- !
957         !                                           ! all cases !
958         !                                           ! --------- !
[9367]959         CALL iom_rstput( kt, nitrst, numrow, 'e3t_b', e3t_b(:,:,:), ldxios = lwxios )
960         CALL iom_rstput( kt, nitrst, numrow, 'e3t_n', e3t_n(:,:,:), ldxios = lwxios )
[4292]961         !                                           ! ----------------------- !
962         IF( ln_vvl_ztilde .OR. ln_vvl_layer ) THEN  ! z_tilde and layer cases !
963            !                                        ! ----------------------- !
[9367]964            CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_b', tilde_e3t_b(:,:,:), ldxios = lwxios)
965            CALL iom_rstput( kt, nitrst, numrow, 'tilde_e3t_n', tilde_e3t_n(:,:,:), ldxios = lwxios)
[4292]966         END IF
967         !                                           ! -------------!   
968         IF( ln_vvl_ztilde ) THEN                    ! z_tilde case !
969            !                                        ! ------------ !
[9367]970            CALL iom_rstput( kt, nitrst, numrow, 'hdiv_lf', hdiv_lf(:,:,:), ldxios = lwxios)
[4292]971         ENDIF
[5836]972         !
[9367]973         IF( lwxios ) CALL iom_swap(      cxios_context          )
[4292]974      ENDIF
[5836]975      !
[4292]976   END SUBROUTINE dom_vvl_rst
977
978
979   SUBROUTINE dom_vvl_ctl
980      !!---------------------------------------------------------------------
981      !!                  ***  ROUTINE dom_vvl_ctl  ***
982      !!               
983      !! ** Purpose :   Control the consistency between namelist options
984      !!                for vertical coordinate
985      !!----------------------------------------------------------------------
[5836]986      INTEGER ::   ioptio, ios
987      !!
[4292]988      NAMELIST/nam_vvl/ ln_vvl_zstar, ln_vvl_ztilde, ln_vvl_layer, ln_vvl_ztilde_as_zstar, &
[5836]989         &              ln_vvl_zstar_at_eqtor      , rn_ahe3     , rn_rst_e3t            , &
990         &              rn_lf_cutoff               , rn_zdef_max , ln_vvl_dbg                ! not yet implemented: ln_vvl_kepe
[4292]991      !!----------------------------------------------------------------------
[5836]992      !
[4294]993      REWIND( numnam_ref )              ! Namelist nam_vvl in reference namelist :
994      READ  ( numnam_ref, nam_vvl, IOSTAT = ios, ERR = 901)
[11536]995901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nam_vvl in reference namelist' )
[4294]996      REWIND( numnam_cfg )              ! Namelist nam_vvl in configuration namelist : Parameters of the run
997      READ  ( numnam_cfg, nam_vvl, IOSTAT = ios, ERR = 902 )
[11536]998902   IF( ios >  0 ) CALL ctl_nam ( ios , 'nam_vvl in configuration namelist' )
[4624]999      IF(lwm) WRITE ( numond, nam_vvl )
[5836]1000      !
[4292]1001      IF(lwp) THEN                    ! Namelist print
1002         WRITE(numout,*)
1003         WRITE(numout,*) 'dom_vvl_ctl : choice/control of the variable vertical coordinate'
1004         WRITE(numout,*) '~~~~~~~~~~~'
[9168]1005         WRITE(numout,*) '   Namelist nam_vvl : chose a vertical coordinate'
1006         WRITE(numout,*) '      zstar                      ln_vvl_zstar   = ', ln_vvl_zstar
1007         WRITE(numout,*) '      ztilde                     ln_vvl_ztilde  = ', ln_vvl_ztilde
1008         WRITE(numout,*) '      layer                      ln_vvl_layer   = ', ln_vvl_layer
1009         WRITE(numout,*) '      ztilde as zstar   ln_vvl_ztilde_as_zstar  = ', ln_vvl_ztilde_as_zstar
[4292]1010         WRITE(numout,*) '      ztilde near the equator    ln_vvl_zstar_at_eqtor  = ', ln_vvl_zstar_at_eqtor
[9168]1011         WRITE(numout,*) '      !'
1012         WRITE(numout,*) '      thickness diffusion coefficient                      rn_ahe3      = ', rn_ahe3
1013         WRITE(numout,*) '      maximum e3t deformation fractional change            rn_zdef_max  = ', rn_zdef_max
[4292]1014         IF( ln_vvl_ztilde_as_zstar ) THEN
[9168]1015            WRITE(numout,*) '      ztilde running in zstar emulation mode (ln_vvl_ztilde_as_zstar=T) '
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.e0'
1019            WRITE(numout,*) '            hard-wired : z-tilde cutoff frequency of low-pass filter (days)'
1020            WRITE(numout,*) '                         rn_lf_cutoff   = 1.0/rdt'
[4292]1021         ELSE
[9168]1022            WRITE(numout,*) '      z-tilde to zstar restoration timescale (days)        rn_rst_e3t   = ', rn_rst_e3t
1023            WRITE(numout,*) '      z-tilde cutoff frequency of low-pass filter (days)   rn_lf_cutoff = ', rn_lf_cutoff
[4292]1024         ENDIF
[9168]1025         WRITE(numout,*) '         debug prints flag                                 ln_vvl_dbg   = ', ln_vvl_dbg
[4292]1026      ENDIF
[5836]1027      !
[4292]1028      ioptio = 0                      ! Parameter control
[5836]1029      IF( ln_vvl_ztilde_as_zstar )   ln_vvl_ztilde = .true.
1030      IF( ln_vvl_zstar           )   ioptio = ioptio + 1
1031      IF( ln_vvl_ztilde          )   ioptio = ioptio + 1
1032      IF( ln_vvl_layer           )   ioptio = ioptio + 1
1033      !
[4292]1034      IF( ioptio /= 1 )   CALL ctl_stop( 'Choose ONE vertical coordinate in namelist nam_vvl' )
[6140]1035      IF( .NOT. ln_vvl_zstar .AND. ln_isf ) CALL ctl_stop( 'Only vvl_zstar has been tested with ice shelf cavity' )
[5836]1036      !
[4292]1037      IF(lwp) THEN                   ! Print the choice
1038         WRITE(numout,*)
[9168]1039         IF( ln_vvl_zstar           ) WRITE(numout,*) '      ==>>>   zstar vertical coordinate is used'
1040         IF( ln_vvl_ztilde          ) WRITE(numout,*) '      ==>>>   ztilde vertical coordinate is used'
1041         IF( ln_vvl_layer           ) WRITE(numout,*) '      ==>>>   layer vertical coordinate is used'
1042         IF( ln_vvl_ztilde_as_zstar ) WRITE(numout,*) '      ==>>>   to emulate a zstar coordinate'
[4292]1043      ENDIF
[5836]1044      !
[4486]1045#if defined key_agrif
[9190]1046      IF( (.NOT.Agrif_Root()).AND.(.NOT.ln_vvl_zstar) )   CALL ctl_stop( 'AGRIF is implemented with zstar coordinate only' )
[4486]1047#endif
[5836]1048      !
[4292]1049   END SUBROUTINE dom_vvl_ctl
1050
[592]1051   !!======================================================================
1052END MODULE domvvl
Note: See TracBrowser for help on using the repository browser.