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 @ 14433

Last change on this file since 14433 was 14433, checked in by smasson, 3 years ago

trunk: merge dev_r14312_MPI_Interface into the trunk, #2598

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